xref: /petsc/src/ts/interface/sensitivity/tssen.c (revision b1e111eb0d8a8c8f79792012adb4e9519947d865)
1814e85d6SHong Zhang #include <petsc/private/tsimpl.h>        /*I "petscts.h"  I*/
2814e85d6SHong Zhang #include <petscdraw.h>
3814e85d6SHong Zhang 
433f72c5dSHong Zhang PetscLogEvent TS_AdjointStep,TS_ForwardStep,TS_JacobianPEval;
5814e85d6SHong Zhang 
6814e85d6SHong Zhang /* ------------------------ Sensitivity Context ---------------------------*/
7814e85d6SHong Zhang 
8814e85d6SHong Zhang /*@C
9b5b4017aSHong Zhang   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.
10814e85d6SHong Zhang 
11814e85d6SHong Zhang   Logically Collective on TS
12814e85d6SHong Zhang 
13814e85d6SHong Zhang   Input Parameters:
14b5b4017aSHong Zhang + ts - TS context obtained from TSCreate()
15b5b4017aSHong Zhang . Amat - JacobianP matrix
16b5b4017aSHong Zhang . func - function
17b5b4017aSHong Zhang - ctx - [optional] user-defined function context
18814e85d6SHong Zhang 
19814e85d6SHong Zhang   Calling sequence of func:
20814e85d6SHong Zhang $ func (TS ts,PetscReal t,Vec y,Mat A,void *ctx);
21814e85d6SHong Zhang +   t - current timestep
22c9ad9de0SHong Zhang .   U - input vector (current ODE solution)
23814e85d6SHong Zhang .   A - output matrix
24814e85d6SHong Zhang -   ctx - [optional] user-defined function context
25814e85d6SHong Zhang 
26814e85d6SHong Zhang   Level: intermediate
27814e85d6SHong Zhang 
28814e85d6SHong Zhang   Notes:
29814e85d6SHong Zhang     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
30814e85d6SHong Zhang 
31814e85d6SHong Zhang .keywords: TS, sensitivity
32cd4cee2dSHong Zhang .seealso: TSGetRHSJacobianP()
33814e85d6SHong Zhang @*/
34814e85d6SHong Zhang PetscErrorCode TSSetRHSJacobianP(TS ts,Mat Amat,PetscErrorCode (*func)(TS,PetscReal,Vec,Mat,void*),void *ctx)
35814e85d6SHong Zhang {
36814e85d6SHong Zhang   PetscErrorCode ierr;
37814e85d6SHong Zhang 
38814e85d6SHong Zhang   PetscFunctionBegin;
39814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
40814e85d6SHong Zhang   PetscValidHeaderSpecific(Amat,MAT_CLASSID,2);
41814e85d6SHong Zhang 
42814e85d6SHong Zhang   ts->rhsjacobianp    = func;
43814e85d6SHong Zhang   ts->rhsjacobianpctx = ctx;
44814e85d6SHong Zhang   if(Amat) {
45814e85d6SHong Zhang     ierr = PetscObjectReference((PetscObject)Amat);CHKERRQ(ierr);
4633f72c5dSHong Zhang     ierr = MatDestroy(&ts->Jacprhs);CHKERRQ(ierr);
4733f72c5dSHong Zhang     ts->Jacprhs = Amat;
48814e85d6SHong Zhang   }
49814e85d6SHong Zhang   PetscFunctionReturn(0);
50814e85d6SHong Zhang }
51814e85d6SHong Zhang 
52814e85d6SHong Zhang /*@C
53cd4cee2dSHong Zhang   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.
54cd4cee2dSHong Zhang 
55cd4cee2dSHong Zhang   Logically Collective on TS
56cd4cee2dSHong Zhang 
57cd4cee2dSHong Zhang   Input Parameters:
58cd4cee2dSHong Zhang . ts - TS context obtained from TSCreate()
59cd4cee2dSHong Zhang 
60cd4cee2dSHong Zhang   Output Parameters:
61cd4cee2dSHong Zhang + Amat - JacobianP matrix
62cd4cee2dSHong Zhang . func - function
63cd4cee2dSHong Zhang - ctx - [optional] user-defined function context
64cd4cee2dSHong Zhang 
65cd4cee2dSHong Zhang   Calling sequence of func:
66cd4cee2dSHong Zhang $ func (TS ts,PetscReal t,Vec y,Mat A,void *ctx);
67cd4cee2dSHong Zhang +   t - current timestep
68cd4cee2dSHong Zhang .   U - input vector (current ODE solution)
69cd4cee2dSHong Zhang .   A - output matrix
70cd4cee2dSHong Zhang -   ctx - [optional] user-defined function context
71cd4cee2dSHong Zhang 
72cd4cee2dSHong Zhang   Level: intermediate
73cd4cee2dSHong Zhang 
74cd4cee2dSHong Zhang   Notes:
75cd4cee2dSHong Zhang     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
76cd4cee2dSHong Zhang 
77cd4cee2dSHong Zhang .keywords: TS, sensitivity
78cd4cee2dSHong Zhang .seealso: TSSetRHSJacobianP()
79cd4cee2dSHong Zhang @*/
80cd4cee2dSHong Zhang PetscErrorCode TSGetRHSJacobianP(TS ts,Mat *Amat,PetscErrorCode (**func)(TS,PetscReal,Vec,Mat,void*),void **ctx)
81cd4cee2dSHong Zhang {
82cd4cee2dSHong Zhang   PetscFunctionBegin;
83cd4cee2dSHong Zhang   if (func) *func = ts->rhsjacobianp;
84cd4cee2dSHong Zhang   if (ctx) *ctx  = ts->rhsjacobianpctx;
85cd4cee2dSHong Zhang   if (Amat) *Amat = ts->Jacprhs;
86cd4cee2dSHong Zhang   PetscFunctionReturn(0);
87cd4cee2dSHong Zhang }
88cd4cee2dSHong Zhang 
89cd4cee2dSHong Zhang /*@C
90814e85d6SHong Zhang   TSComputeRHSJacobianP - Runs the user-defined JacobianP function.
91814e85d6SHong Zhang 
92814e85d6SHong Zhang   Collective on TS
93814e85d6SHong Zhang 
94814e85d6SHong Zhang   Input Parameters:
95814e85d6SHong Zhang . ts   - The TS context obtained from TSCreate()
96814e85d6SHong Zhang 
97814e85d6SHong Zhang   Level: developer
98814e85d6SHong Zhang 
99814e85d6SHong Zhang .keywords: TS, sensitivity
100814e85d6SHong Zhang .seealso: TSSetRHSJacobianP()
101814e85d6SHong Zhang @*/
102c9ad9de0SHong Zhang PetscErrorCode TSComputeRHSJacobianP(TS ts,PetscReal t,Vec U,Mat Amat)
103814e85d6SHong Zhang {
104814e85d6SHong Zhang   PetscErrorCode ierr;
105814e85d6SHong Zhang 
106814e85d6SHong Zhang   PetscFunctionBegin;
10733f72c5dSHong Zhang   if (!Amat) PetscFunctionReturn(0);
108814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
109c9ad9de0SHong Zhang   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
110814e85d6SHong Zhang 
111814e85d6SHong Zhang   PetscStackPush("TS user JacobianP function for sensitivity analysis");
112c9ad9de0SHong Zhang   ierr = (*ts->rhsjacobianp)(ts,t,U,Amat,ts->rhsjacobianpctx);CHKERRQ(ierr);
113814e85d6SHong Zhang   PetscStackPop;
114814e85d6SHong Zhang   PetscFunctionReturn(0);
115814e85d6SHong Zhang }
116814e85d6SHong Zhang 
117814e85d6SHong Zhang /*@C
11833f72c5dSHong Zhang   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.
11933f72c5dSHong Zhang 
12033f72c5dSHong Zhang   Logically Collective on TS
12133f72c5dSHong Zhang 
12233f72c5dSHong Zhang   Input Parameters:
12333f72c5dSHong Zhang + ts - TS context obtained from TSCreate()
12433f72c5dSHong Zhang . Amat - JacobianP matrix
12533f72c5dSHong Zhang . func - function
12633f72c5dSHong Zhang - ctx - [optional] user-defined function context
12733f72c5dSHong Zhang 
12833f72c5dSHong Zhang   Calling sequence of func:
12933f72c5dSHong Zhang $ func (TS ts,PetscReal t,Vec y,Mat A,void *ctx);
13033f72c5dSHong Zhang +   t - current timestep
13133f72c5dSHong Zhang .   U - input vector (current ODE solution)
13233f72c5dSHong Zhang .   Udot - time derivative of state vector
13333f72c5dSHong Zhang .   shift - shift to apply, see note below
13433f72c5dSHong Zhang .   A - output matrix
13533f72c5dSHong Zhang -   ctx - [optional] user-defined function context
13633f72c5dSHong Zhang 
13733f72c5dSHong Zhang   Level: intermediate
13833f72c5dSHong Zhang 
13933f72c5dSHong Zhang   Notes:
14033f72c5dSHong Zhang     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
14133f72c5dSHong Zhang 
14233f72c5dSHong Zhang .keywords: TS, sensitivity
14333f72c5dSHong Zhang .seealso:
14433f72c5dSHong Zhang @*/
14533f72c5dSHong Zhang PetscErrorCode TSSetIJacobianP(TS ts,Mat Amat,PetscErrorCode (*func)(TS,PetscReal,Vec,Vec,PetscReal,Mat,void*),void *ctx)
14633f72c5dSHong Zhang {
14733f72c5dSHong Zhang   PetscErrorCode ierr;
14833f72c5dSHong Zhang 
14933f72c5dSHong Zhang   PetscFunctionBegin;
15033f72c5dSHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
15133f72c5dSHong Zhang   PetscValidHeaderSpecific(Amat,MAT_CLASSID,2);
15233f72c5dSHong Zhang 
15333f72c5dSHong Zhang   ts->ijacobianp    = func;
15433f72c5dSHong Zhang   ts->ijacobianpctx = ctx;
15533f72c5dSHong Zhang   if(Amat) {
15633f72c5dSHong Zhang     ierr = PetscObjectReference((PetscObject)Amat);CHKERRQ(ierr);
15733f72c5dSHong Zhang     ierr = MatDestroy(&ts->Jacp);CHKERRQ(ierr);
15833f72c5dSHong Zhang     ts->Jacp = Amat;
15933f72c5dSHong Zhang   }
16033f72c5dSHong Zhang   PetscFunctionReturn(0);
16133f72c5dSHong Zhang }
16233f72c5dSHong Zhang 
16333f72c5dSHong Zhang /*@C
16433f72c5dSHong Zhang   TSComputeIJacobianP - Runs the user-defined IJacobianP function.
16533f72c5dSHong Zhang 
16633f72c5dSHong Zhang   Collective on TS
16733f72c5dSHong Zhang 
16833f72c5dSHong Zhang   Input Parameters:
16933f72c5dSHong Zhang + ts - the TS context
17033f72c5dSHong Zhang . t - current timestep
17133f72c5dSHong Zhang . U - state vector
17233f72c5dSHong Zhang . Udot - time derivative of state vector
17333f72c5dSHong Zhang . shift - shift to apply, see note below
17433f72c5dSHong Zhang - imex - flag indicates if the method is IMEX so that the RHSJacobian should be kept separate
17533f72c5dSHong Zhang 
17633f72c5dSHong Zhang   Output Parameters:
17733f72c5dSHong Zhang . A - Jacobian matrix
17833f72c5dSHong Zhang 
17933f72c5dSHong Zhang   Level: developer
18033f72c5dSHong Zhang 
18133f72c5dSHong Zhang .keywords: TS, sensitivity
18233f72c5dSHong Zhang .seealso: TSSetIJacobianP()
18333f72c5dSHong Zhang @*/
18433f72c5dSHong Zhang PetscErrorCode TSComputeIJacobianP(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal shift,Mat Amat,PetscBool imex)
18533f72c5dSHong Zhang {
18633f72c5dSHong Zhang   PetscErrorCode ierr;
18733f72c5dSHong Zhang 
18833f72c5dSHong Zhang   PetscFunctionBegin;
18933f72c5dSHong Zhang   if (!Amat) PetscFunctionReturn(0);
19033f72c5dSHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
19133f72c5dSHong Zhang   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
19233f72c5dSHong Zhang   PetscValidHeaderSpecific(Udot,VEC_CLASSID,4);
19333f72c5dSHong Zhang 
19433f72c5dSHong Zhang   ierr = PetscLogEventBegin(TS_JacobianPEval,ts,U,Amat,0);CHKERRQ(ierr);
19533f72c5dSHong Zhang   if (ts->ijacobianp) {
19633f72c5dSHong Zhang     PetscStackPush("TS user JacobianP function for sensitivity analysis");
19733f72c5dSHong Zhang     ierr = (*ts->ijacobianp)(ts,t,U,Udot,shift,Amat,ts->ijacobianpctx);CHKERRQ(ierr);
19833f72c5dSHong Zhang     PetscStackPop;
19933f72c5dSHong Zhang   }
20033f72c5dSHong Zhang   if (imex) {
20133f72c5dSHong Zhang     if (!ts->ijacobianp) {  /* system was written as Udot = G(t,U) */
20233f72c5dSHong Zhang       PetscBool assembled;
20333f72c5dSHong Zhang       ierr = MatZeroEntries(Amat);CHKERRQ(ierr);
20433f72c5dSHong Zhang       ierr = MatAssembled(Amat,&assembled);CHKERRQ(ierr);
20533f72c5dSHong Zhang       if (!assembled) {
20633f72c5dSHong Zhang         ierr = MatAssemblyBegin(Amat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
20733f72c5dSHong Zhang         ierr = MatAssemblyEnd(Amat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
20833f72c5dSHong Zhang       }
20933f72c5dSHong Zhang     }
21033f72c5dSHong Zhang   } else {
21133f72c5dSHong Zhang     if (ts->rhsjacobianp) {
21233f72c5dSHong Zhang       ierr = TSComputeRHSJacobianP(ts,t,U,ts->Jacprhs);CHKERRQ(ierr);
21333f72c5dSHong Zhang     }
21433f72c5dSHong Zhang     if (ts->Jacprhs == Amat) { /* No IJacobian, so we only have the RHS matrix */
21533f72c5dSHong Zhang       ierr = MatScale(Amat,-1);CHKERRQ(ierr);
21633f72c5dSHong Zhang     } else if (ts->Jacprhs) { /* Both IJacobian and RHSJacobian */
21733f72c5dSHong Zhang       MatStructure axpy = DIFFERENT_NONZERO_PATTERN;
21833f72c5dSHong Zhang       if (!ts->ijacobianp) { /* No IJacobianp provided, but we have a separate RHS matrix */
21933f72c5dSHong Zhang         ierr = MatZeroEntries(Amat);CHKERRQ(ierr);
22033f72c5dSHong Zhang       }
22133f72c5dSHong Zhang       ierr = MatAXPY(Amat,-1,ts->Jacprhs,axpy);CHKERRQ(ierr);
22233f72c5dSHong Zhang     }
22333f72c5dSHong Zhang   }
22433f72c5dSHong Zhang   ierr = PetscLogEventEnd(TS_JacobianPEval,ts,U,Amat,0);CHKERRQ(ierr);
22533f72c5dSHong Zhang   PetscFunctionReturn(0);
22633f72c5dSHong Zhang }
22733f72c5dSHong Zhang 
22833f72c5dSHong Zhang /*@C
229814e85d6SHong Zhang     TSSetCostIntegrand - Sets the routine for evaluating the integral term in one or more cost functions
230814e85d6SHong Zhang 
231814e85d6SHong Zhang     Logically Collective on TS
232814e85d6SHong Zhang 
233814e85d6SHong Zhang     Input Parameters:
234814e85d6SHong Zhang +   ts - the TS context obtained from TSCreate()
235814e85d6SHong Zhang .   numcost - number of gradients to be computed, this is the number of cost functions
236814e85d6SHong Zhang .   costintegral - vector that stores the integral values
237814e85d6SHong Zhang .   rf - routine for evaluating the integrand function
238c9ad9de0SHong Zhang .   drduf - function that computes the gradients of the r's with respect to u
239814e85d6SHong Zhang .   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)
240814e85d6SHong Zhang .   fwd - flag indicating whether to evaluate cost integral in the forward run or the adjoint run
241814e85d6SHong Zhang -   ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
242814e85d6SHong Zhang 
243814e85d6SHong Zhang     Calling sequence of rf:
244c9ad9de0SHong Zhang $   PetscErrorCode rf(TS ts,PetscReal t,Vec U,Vec F,void *ctx);
245814e85d6SHong Zhang 
246c9ad9de0SHong Zhang     Calling sequence of drduf:
247c9ad9de0SHong Zhang $   PetscErroCode drduf(TS ts,PetscReal t,Vec U,Vec *dRdU,void *ctx);
248814e85d6SHong Zhang 
249814e85d6SHong Zhang     Calling sequence of drdpf:
250c9ad9de0SHong Zhang $   PetscErroCode drdpf(TS ts,PetscReal t,Vec U,Vec *dRdP,void *ctx);
251814e85d6SHong Zhang 
252cd4cee2dSHong Zhang     Level: deprecated
253814e85d6SHong Zhang 
254814e85d6SHong Zhang     Notes:
255814e85d6SHong Zhang     For optimization there is usually a single cost function (numcost = 1). For sensitivities there may be multiple cost functions
256814e85d6SHong Zhang 
257814e85d6SHong Zhang .keywords: TS, sensitivity analysis, timestep, set, quadrature, function
258814e85d6SHong Zhang 
259814e85d6SHong Zhang .seealso: TSSetRHSJacobianP(), TSGetCostGradients(), TSSetCostGradients()
260814e85d6SHong Zhang @*/
261814e85d6SHong Zhang PetscErrorCode TSSetCostIntegrand(TS ts,PetscInt numcost,Vec costintegral,PetscErrorCode (*rf)(TS,PetscReal,Vec,Vec,void*),
262c9ad9de0SHong Zhang                                                           PetscErrorCode (*drduf)(TS,PetscReal,Vec,Vec*,void*),
263814e85d6SHong Zhang                                                           PetscErrorCode (*drdpf)(TS,PetscReal,Vec,Vec*,void*),
264814e85d6SHong Zhang                                                           PetscBool fwd,void *ctx)
265814e85d6SHong Zhang {
266814e85d6SHong Zhang   PetscErrorCode ierr;
267814e85d6SHong Zhang 
268814e85d6SHong Zhang   PetscFunctionBegin;
269814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
270814e85d6SHong Zhang   if (costintegral) PetscValidHeaderSpecific(costintegral,VEC_CLASSID,3);
271814e85d6SHong Zhang   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()");
272814e85d6SHong Zhang   if (!ts->numcost) ts->numcost=numcost;
273814e85d6SHong Zhang 
274814e85d6SHong Zhang   if (costintegral) {
275814e85d6SHong Zhang     ierr = PetscObjectReference((PetscObject)costintegral);CHKERRQ(ierr);
276814e85d6SHong Zhang     ierr = VecDestroy(&ts->vec_costintegral);CHKERRQ(ierr);
277814e85d6SHong Zhang     ts->vec_costintegral = costintegral;
278814e85d6SHong Zhang   } else {
279814e85d6SHong Zhang     if (!ts->vec_costintegral) { /* Create a seq vec if user does not provide one */
280814e85d6SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,numcost,&ts->vec_costintegral);CHKERRQ(ierr);
281814e85d6SHong Zhang     } else {
282814e85d6SHong Zhang       ierr = VecSet(ts->vec_costintegral,0.0);CHKERRQ(ierr);
283814e85d6SHong Zhang     }
284814e85d6SHong Zhang   }
285814e85d6SHong Zhang   if (!ts->vec_costintegrand) {
286814e85d6SHong Zhang     ierr = VecDuplicate(ts->vec_costintegral,&ts->vec_costintegrand);CHKERRQ(ierr);
287814e85d6SHong Zhang   } else {
288814e85d6SHong Zhang     ierr = VecSet(ts->vec_costintegrand,0.0);CHKERRQ(ierr);
289814e85d6SHong Zhang   }
290814e85d6SHong Zhang   ts->costintegralfwd  = fwd; /* Evaluate the cost integral in forward run if fwd is true */
291814e85d6SHong Zhang   ts->costintegrand    = rf;
292814e85d6SHong Zhang   ts->costintegrandctx = ctx;
293c9ad9de0SHong Zhang   ts->drdufunction     = drduf;
294814e85d6SHong Zhang   ts->drdpfunction     = drdpf;
295814e85d6SHong Zhang   PetscFunctionReturn(0);
296814e85d6SHong Zhang }
297814e85d6SHong Zhang 
298b5b4017aSHong Zhang /*@C
299814e85d6SHong Zhang    TSGetCostIntegral - Returns the values of the integral term in the cost functions.
300814e85d6SHong Zhang    It is valid to call the routine after a backward run.
301814e85d6SHong Zhang 
302814e85d6SHong Zhang    Not Collective
303814e85d6SHong Zhang 
304814e85d6SHong Zhang    Input Parameter:
305814e85d6SHong Zhang .  ts - the TS context obtained from TSCreate()
306814e85d6SHong Zhang 
307814e85d6SHong Zhang    Output Parameter:
308814e85d6SHong Zhang .  v - the vector containing the integrals for each cost function
309814e85d6SHong Zhang 
310814e85d6SHong Zhang    Level: intermediate
311814e85d6SHong Zhang 
312814e85d6SHong Zhang .seealso: TSSetCostIntegrand()
313814e85d6SHong Zhang 
314814e85d6SHong Zhang .keywords: TS, sensitivity analysis
315814e85d6SHong Zhang @*/
316814e85d6SHong Zhang PetscErrorCode  TSGetCostIntegral(TS ts,Vec *v)
317814e85d6SHong Zhang {
318cd4cee2dSHong Zhang   TS             quadts;
319cd4cee2dSHong Zhang   PetscErrorCode ierr;
320cd4cee2dSHong Zhang 
321814e85d6SHong Zhang   PetscFunctionBegin;
322814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
323814e85d6SHong Zhang   PetscValidPointer(v,2);
324cd4cee2dSHong Zhang   ierr = TSGetQuadratureTS(ts,NULL,&quadts);CHKERRQ(ierr);
325cd4cee2dSHong Zhang   *v = quadts->vec_sol;
326814e85d6SHong Zhang   PetscFunctionReturn(0);
327814e85d6SHong Zhang }
328814e85d6SHong Zhang 
329b5b4017aSHong Zhang /*@C
330814e85d6SHong Zhang    TSComputeCostIntegrand - Evaluates the integral function in the cost functions.
331814e85d6SHong Zhang 
332814e85d6SHong Zhang    Input Parameters:
333814e85d6SHong Zhang +  ts - the TS context
334814e85d6SHong Zhang .  t - current time
335c9ad9de0SHong Zhang -  U - state vector, i.e. current solution
336814e85d6SHong Zhang 
337814e85d6SHong Zhang    Output Parameter:
338c9ad9de0SHong Zhang .  Q - vector of size numcost to hold the outputs
339814e85d6SHong Zhang 
340369cf35fSHong Zhang    Notes:
341814e85d6SHong Zhang    Most users should not need to explicitly call this routine, as it
342814e85d6SHong Zhang    is used internally within the sensitivity analysis context.
343814e85d6SHong Zhang 
344cd4cee2dSHong Zhang    Level: deprecated
345814e85d6SHong Zhang 
346814e85d6SHong Zhang .keywords: TS, compute
347814e85d6SHong Zhang 
348814e85d6SHong Zhang .seealso: TSSetCostIntegrand()
349814e85d6SHong Zhang @*/
350c9ad9de0SHong Zhang PetscErrorCode TSComputeCostIntegrand(TS ts,PetscReal t,Vec U,Vec Q)
351814e85d6SHong Zhang {
352814e85d6SHong Zhang   PetscErrorCode ierr;
353814e85d6SHong Zhang 
354814e85d6SHong Zhang   PetscFunctionBegin;
355814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
356c9ad9de0SHong Zhang   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
357c9ad9de0SHong Zhang   PetscValidHeaderSpecific(Q,VEC_CLASSID,4);
358814e85d6SHong Zhang 
359c9ad9de0SHong Zhang   ierr = PetscLogEventBegin(TS_FunctionEval,ts,U,Q,0);CHKERRQ(ierr);
360814e85d6SHong Zhang   if (ts->costintegrand) {
361814e85d6SHong Zhang     PetscStackPush("TS user integrand in the cost function");
362c9ad9de0SHong Zhang     ierr = (*ts->costintegrand)(ts,t,U,Q,ts->costintegrandctx);CHKERRQ(ierr);
363814e85d6SHong Zhang     PetscStackPop;
364814e85d6SHong Zhang   } else {
365c9ad9de0SHong Zhang     ierr = VecZeroEntries(Q);CHKERRQ(ierr);
366814e85d6SHong Zhang   }
367814e85d6SHong Zhang 
368c9ad9de0SHong Zhang   ierr = PetscLogEventEnd(TS_FunctionEval,ts,U,Q,0);CHKERRQ(ierr);
369814e85d6SHong Zhang   PetscFunctionReturn(0);
370814e85d6SHong Zhang }
371814e85d6SHong Zhang 
372b5b4017aSHong Zhang /*@C
373c9ad9de0SHong Zhang   TSComputeDRDUFunction - Runs the user-defined DRDU function.
374814e85d6SHong Zhang 
375814e85d6SHong Zhang   Collective on TS
376814e85d6SHong Zhang 
377814e85d6SHong Zhang   Input Parameters:
378c9ad9de0SHong Zhang + ts - the TS context obtained from TSCreate()
379c9ad9de0SHong Zhang . t - current time
380c9ad9de0SHong Zhang - U - stata vector
381c9ad9de0SHong Zhang 
382c9ad9de0SHong Zhang   Output Parameters:
383b5b4017aSHong Zhang . DRDU - vector array to hold the outputs
384814e85d6SHong Zhang 
385814e85d6SHong Zhang   Notes:
386c9ad9de0SHong Zhang   TSComputeDRDUFunction() is typically used for sensitivity implementation,
387814e85d6SHong Zhang   so most users would not generally call this routine themselves.
388814e85d6SHong Zhang 
389814e85d6SHong Zhang   Level: developer
390814e85d6SHong Zhang 
391814e85d6SHong Zhang .keywords: TS, sensitivity
392814e85d6SHong Zhang .seealso: TSSetCostIntegrand()
393814e85d6SHong Zhang @*/
394c9ad9de0SHong Zhang PetscErrorCode TSComputeDRDUFunction(TS ts,PetscReal t,Vec U,Vec *DRDU)
395814e85d6SHong Zhang {
396814e85d6SHong Zhang   PetscErrorCode ierr;
397814e85d6SHong Zhang 
398814e85d6SHong Zhang   PetscFunctionBegin;
39933f72c5dSHong Zhang   if (!DRDU) PetscFunctionReturn(0);
400814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
401c9ad9de0SHong Zhang   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
402814e85d6SHong Zhang 
403c9ad9de0SHong Zhang   PetscStackPush("TS user DRDU function for sensitivity analysis");
404c9ad9de0SHong Zhang   ierr = (*ts->drdufunction)(ts,t,U,DRDU,ts->costintegrandctx); CHKERRQ(ierr);
405814e85d6SHong Zhang   PetscStackPop;
406814e85d6SHong Zhang   PetscFunctionReturn(0);
407814e85d6SHong Zhang }
408814e85d6SHong Zhang 
409b5b4017aSHong Zhang /*@C
410814e85d6SHong Zhang   TSComputeDRDPFunction - Runs the user-defined DRDP function.
411814e85d6SHong Zhang 
412814e85d6SHong Zhang   Collective on TS
413814e85d6SHong Zhang 
414814e85d6SHong Zhang   Input Parameters:
415c9ad9de0SHong Zhang + ts - the TS context obtained from TSCreate()
416c9ad9de0SHong Zhang . t - current time
417c9ad9de0SHong Zhang - U - stata vector
418c9ad9de0SHong Zhang 
419c9ad9de0SHong Zhang   Output Parameters:
420b5b4017aSHong Zhang . DRDP - vector array to hold the outputs
421814e85d6SHong Zhang 
422814e85d6SHong Zhang   Notes:
423c9ad9de0SHong Zhang   TSComputeDRDPFunction() is typically used for sensitivity implementation,
424814e85d6SHong Zhang   so most users would not generally call this routine themselves.
425814e85d6SHong Zhang 
426814e85d6SHong Zhang   Level: developer
427814e85d6SHong Zhang 
428814e85d6SHong Zhang .keywords: TS, sensitivity
429814e85d6SHong Zhang .seealso: TSSetCostIntegrand()
430814e85d6SHong Zhang @*/
431c9ad9de0SHong Zhang PetscErrorCode TSComputeDRDPFunction(TS ts,PetscReal t,Vec U,Vec *DRDP)
432814e85d6SHong Zhang {
433814e85d6SHong Zhang   PetscErrorCode ierr;
434814e85d6SHong Zhang 
435814e85d6SHong Zhang   PetscFunctionBegin;
43633f72c5dSHong Zhang   if (!DRDP) PetscFunctionReturn(0);
437814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
438c9ad9de0SHong Zhang   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
439814e85d6SHong Zhang 
440814e85d6SHong Zhang   PetscStackPush("TS user DRDP function for sensitivity analysis");
441c9ad9de0SHong Zhang   ierr = (*ts->drdpfunction)(ts,t,U,DRDP,ts->costintegrandctx); CHKERRQ(ierr);
442814e85d6SHong Zhang   PetscStackPop;
443814e85d6SHong Zhang   PetscFunctionReturn(0);
444814e85d6SHong Zhang }
445814e85d6SHong Zhang 
446b5b4017aSHong Zhang /*@C
447b5b4017aSHong Zhang   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.
448b5b4017aSHong Zhang 
449b5b4017aSHong Zhang   Logically Collective on TS
450b5b4017aSHong Zhang 
451b5b4017aSHong Zhang   Input Parameters:
452b5b4017aSHong Zhang + ts - TS context obtained from TSCreate()
453b5b4017aSHong Zhang . ihp1 - an array of vectors storing the result of vector-Hessian-vector product for F_UU
454b5b4017aSHong Zhang . hessianproductfunc1 - vector-Hessian-vector product function for F_UU
455b5b4017aSHong Zhang . ihp2 - an array of vectors storing the result of vector-Hessian-vector product for F_UP
456b5b4017aSHong Zhang . hessianproductfunc2 - vector-Hessian-vector product function for F_UP
457b5b4017aSHong Zhang . ihp3 - an array of vectors storing the result of vector-Hessian-vector product for F_PU
458b5b4017aSHong Zhang . hessianproductfunc3 - vector-Hessian-vector product function for F_PU
459b5b4017aSHong Zhang . ihp4 - an array of vectors storing the result of vector-Hessian-vector product for F_PP
460b5b4017aSHong Zhang . hessianproductfunc4 - vector-Hessian-vector product function for F_PP
461b5b4017aSHong Zhang 
462b5b4017aSHong Zhang   Calling sequence of ihessianproductfunc:
463369cf35fSHong Zhang $ ihessianproductfunc (TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx);
464b5b4017aSHong Zhang +   t - current timestep
465b5b4017aSHong Zhang .   U - input vector (current ODE solution)
466369cf35fSHong Zhang .   Vl - an array of input vectors to be left-multiplied with the Hessian
467b5b4017aSHong Zhang .   Vr - input vector to be right-multiplied with the Hessian
468369cf35fSHong Zhang .   VHV - an array of output vectors for vector-Hessian-vector product
469b5b4017aSHong Zhang -   ctx - [optional] user-defined function context
470b5b4017aSHong Zhang 
471b5b4017aSHong Zhang   Level: intermediate
472b5b4017aSHong Zhang 
473369cf35fSHong Zhang   Notes:
474369cf35fSHong Zhang   The first Hessian function and the working array are required.
475369cf35fSHong Zhang   As an example to implement the callback functions, the second callback function calculates the vector-Hessian-vector product
476369cf35fSHong Zhang   $ Vl_n^T*F_UP*Vr
477369cf35fSHong Zhang   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.
478369cf35fSHong Zhang   Each entry of F_UP corresponds to the derivative
479369cf35fSHong Zhang   $ F_UP[i][j][k] = \frac{\partial^2 F[i]}{\partial U[j] \partial P[k]}.
480369cf35fSHong Zhang   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
481369cf35fSHong Zhang   $ VHV_n[j] = \sum_i \sum_k {Vl_n[i] * F_UP[i][j][k] * Vr[k]}
482369cf35fSHong Zhang   If the cost function is a scalar, there will be only one vector in Vl and VHV.
483b5b4017aSHong Zhang 
484b5b4017aSHong Zhang .keywords: TS, sensitivity
485b5b4017aSHong Zhang 
486b5b4017aSHong Zhang .seealso:
487b5b4017aSHong Zhang @*/
488b5b4017aSHong Zhang PetscErrorCode TSSetIHessianProduct(TS ts,Vec *ihp1,PetscErrorCode (*ihessianproductfunc1)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*),
489b5b4017aSHong Zhang                                           Vec *ihp2,PetscErrorCode (*ihessianproductfunc2)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*),
490b5b4017aSHong Zhang                                           Vec *ihp3,PetscErrorCode (*ihessianproductfunc3)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*),
491b5b4017aSHong Zhang                                           Vec *ihp4,PetscErrorCode (*ihessianproductfunc4)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*),
492b5b4017aSHong Zhang                                     void *ctx)
493b5b4017aSHong Zhang {
494b5b4017aSHong Zhang   PetscFunctionBegin;
495b5b4017aSHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
496b5b4017aSHong Zhang   PetscValidPointer(ihp1,2);
497b5b4017aSHong Zhang 
498b5b4017aSHong Zhang   ts->ihessianproductctx = ctx;
499b5b4017aSHong Zhang   if (ihp1) ts->vecs_fuu = ihp1;
500b5b4017aSHong Zhang   if (ihp2) ts->vecs_fup = ihp2;
501b5b4017aSHong Zhang   if (ihp3) ts->vecs_fpu = ihp3;
502b5b4017aSHong Zhang   if (ihp4) ts->vecs_fpp = ihp4;
503b5b4017aSHong Zhang   ts->ihessianproduct_fuu = ihessianproductfunc1;
504b5b4017aSHong Zhang   ts->ihessianproduct_fup = ihessianproductfunc2;
505b5b4017aSHong Zhang   ts->ihessianproduct_fpu = ihessianproductfunc3;
506b5b4017aSHong Zhang   ts->ihessianproduct_fpp = ihessianproductfunc4;
507b5b4017aSHong Zhang   PetscFunctionReturn(0);
508b5b4017aSHong Zhang }
509b5b4017aSHong Zhang 
510b5b4017aSHong Zhang /*@C
511*b1e111ebSHong Zhang   TSComputeIHessianProductFunctionUU - Runs the user-defined vector-Hessian-vector product function for Fuu.
512b5b4017aSHong Zhang 
513b5b4017aSHong Zhang   Collective on TS
514b5b4017aSHong Zhang 
515b5b4017aSHong Zhang   Input Parameters:
516b5b4017aSHong Zhang . ts   - The TS context obtained from TSCreate()
517b5b4017aSHong Zhang 
518b5b4017aSHong Zhang   Notes:
519*b1e111ebSHong Zhang   TSComputeIHessianProductFunctionUU() is typically used for sensitivity implementation,
520b5b4017aSHong Zhang   so most users would not generally call this routine themselves.
521b5b4017aSHong Zhang 
522b5b4017aSHong Zhang   Level: developer
523b5b4017aSHong Zhang 
524b5b4017aSHong Zhang .keywords: TS, sensitivity
525b5b4017aSHong Zhang 
526b5b4017aSHong Zhang .seealso: TSSetIHessianProduct()
527b5b4017aSHong Zhang @*/
528*b1e111ebSHong Zhang PetscErrorCode TSComputeIHessianProductFunctionUU(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV)
529b5b4017aSHong Zhang {
530b5b4017aSHong Zhang   PetscErrorCode ierr;
531b5b4017aSHong Zhang 
532b5b4017aSHong Zhang   PetscFunctionBegin;
53333f72c5dSHong Zhang   if (!VHV) PetscFunctionReturn(0);
534b5b4017aSHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
535b5b4017aSHong Zhang   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
536b5b4017aSHong Zhang 
53767633408SHong Zhang   if (ts->ihessianproduct_fuu) {
538b5b4017aSHong Zhang     PetscStackPush("TS user IHessianProduct function 1 for sensitivity analysis");
539b5b4017aSHong Zhang     ierr = (*ts->ihessianproduct_fuu)(ts,t,U,Vl,Vr,VHV,ts->ihessianproductctx);CHKERRQ(ierr);
540b5b4017aSHong Zhang     PetscStackPop;
54167633408SHong Zhang   }
54267633408SHong Zhang   /* does not consider IMEX for now, so either IHessian or RHSHessian will be calculated, using the same output VHV */
54367633408SHong Zhang   if (ts->rhshessianproduct_guu) {
54467633408SHong Zhang     PetscInt nadj;
545*b1e111ebSHong Zhang     ierr = TSComputeRHSHessianProductFunctionUU(ts,t,U,Vl,Vr,VHV);CHKERRQ(ierr);
54667633408SHong Zhang     for (nadj=0; nadj<ts->numcost; nadj++) {
54767633408SHong Zhang       ierr = VecScale(VHV[nadj],-1);CHKERRQ(ierr);
54867633408SHong Zhang     }
54967633408SHong Zhang   }
550b5b4017aSHong Zhang   PetscFunctionReturn(0);
551b5b4017aSHong Zhang }
552b5b4017aSHong Zhang 
553b5b4017aSHong Zhang /*@C
554*b1e111ebSHong Zhang   TSComputeIHessianProductFunctionUP - Runs the user-defined vector-Hessian-vector product function for Fup.
555b5b4017aSHong Zhang 
556b5b4017aSHong Zhang   Collective on TS
557b5b4017aSHong Zhang 
558b5b4017aSHong Zhang   Input Parameters:
559b5b4017aSHong Zhang . ts   - The TS context obtained from TSCreate()
560b5b4017aSHong Zhang 
561b5b4017aSHong Zhang   Notes:
562*b1e111ebSHong Zhang   TSComputeIHessianProductFunctionUP() is typically used for sensitivity implementation,
563b5b4017aSHong Zhang   so most users would not generally call this routine themselves.
564b5b4017aSHong Zhang 
565b5b4017aSHong Zhang   Level: developer
566b5b4017aSHong Zhang 
567b5b4017aSHong Zhang .keywords: TS, sensitivity
568b5b4017aSHong Zhang 
569b5b4017aSHong Zhang .seealso: TSSetIHessianProduct()
570b5b4017aSHong Zhang @*/
571*b1e111ebSHong Zhang PetscErrorCode TSComputeIHessianProductFunctionUP(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV)
572b5b4017aSHong Zhang {
573b5b4017aSHong Zhang   PetscErrorCode ierr;
574b5b4017aSHong Zhang 
575b5b4017aSHong Zhang   PetscFunctionBegin;
57633f72c5dSHong Zhang   if (!VHV) PetscFunctionReturn(0);
577b5b4017aSHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
578b5b4017aSHong Zhang   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
579b5b4017aSHong Zhang 
58067633408SHong Zhang   if (ts->ihessianproduct_fup) {
581b5b4017aSHong Zhang     PetscStackPush("TS user IHessianProduct function 2 for sensitivity analysis");
582b5b4017aSHong Zhang     ierr = (*ts->ihessianproduct_fup)(ts,t,U,Vl,Vr,VHV,ts->ihessianproductctx);CHKERRQ(ierr);
583b5b4017aSHong Zhang     PetscStackPop;
58467633408SHong Zhang   }
58567633408SHong Zhang   /* does not consider IMEX for now, so either IHessian or RHSHessian will be calculated, using the same output VHV */
58667633408SHong Zhang   if (ts->rhshessianproduct_gup) {
58767633408SHong Zhang     PetscInt nadj;
588*b1e111ebSHong Zhang     ierr = TSComputeRHSHessianProductFunctionUP(ts,t,U,Vl,Vr,VHV);CHKERRQ(ierr);
58967633408SHong Zhang     for (nadj=0; nadj<ts->numcost; nadj++) {
59067633408SHong Zhang       ierr = VecScale(VHV[nadj],-1);CHKERRQ(ierr);
59167633408SHong Zhang     }
59267633408SHong Zhang   }
593b5b4017aSHong Zhang   PetscFunctionReturn(0);
594b5b4017aSHong Zhang }
595b5b4017aSHong Zhang 
596b5b4017aSHong Zhang /*@C
597*b1e111ebSHong Zhang   TSComputeIHessianProductFunctionPU - Runs the user-defined vector-Hessian-vector product function for Fpu.
598b5b4017aSHong Zhang 
599b5b4017aSHong Zhang   Collective on TS
600b5b4017aSHong Zhang 
601b5b4017aSHong Zhang   Input Parameters:
602b5b4017aSHong Zhang . ts   - The TS context obtained from TSCreate()
603b5b4017aSHong Zhang 
604b5b4017aSHong Zhang   Notes:
605*b1e111ebSHong Zhang   TSComputeIHessianProductFunctionPU() is typically used for sensitivity implementation,
606b5b4017aSHong Zhang   so most users would not generally call this routine themselves.
607b5b4017aSHong Zhang 
608b5b4017aSHong Zhang   Level: developer
609b5b4017aSHong Zhang 
610b5b4017aSHong Zhang .keywords: TS, sensitivity
611b5b4017aSHong Zhang 
612b5b4017aSHong Zhang .seealso: TSSetIHessianProduct()
613b5b4017aSHong Zhang @*/
614*b1e111ebSHong Zhang PetscErrorCode TSComputeIHessianProductFunctionPU(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV)
615b5b4017aSHong Zhang {
616b5b4017aSHong Zhang   PetscErrorCode ierr;
617b5b4017aSHong Zhang 
618b5b4017aSHong Zhang   PetscFunctionBegin;
61933f72c5dSHong Zhang   if (!VHV) PetscFunctionReturn(0);
620b5b4017aSHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
621b5b4017aSHong Zhang   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
622b5b4017aSHong Zhang 
62367633408SHong Zhang   if (ts->ihessianproduct_fpu) {
624b5b4017aSHong Zhang     PetscStackPush("TS user IHessianProduct function 3 for sensitivity analysis");
625b5b4017aSHong Zhang     ierr = (*ts->ihessianproduct_fpu)(ts,t,U,Vl,Vr,VHV,ts->ihessianproductctx);CHKERRQ(ierr);
626b5b4017aSHong Zhang     PetscStackPop;
62767633408SHong Zhang   }
62867633408SHong Zhang   /* does not consider IMEX for now, so either IHessian or RHSHessian will be calculated, using the same output VHV */
62967633408SHong Zhang   if (ts->rhshessianproduct_gpu) {
63067633408SHong Zhang     PetscInt nadj;
631*b1e111ebSHong Zhang     ierr = TSComputeRHSHessianProductFunctionPU(ts,t,U,Vl,Vr,VHV);CHKERRQ(ierr);
63267633408SHong Zhang     for (nadj=0; nadj<ts->numcost; nadj++) {
63367633408SHong Zhang       ierr = VecScale(VHV[nadj],-1);CHKERRQ(ierr);
63467633408SHong Zhang     }
63567633408SHong Zhang   }
636b5b4017aSHong Zhang   PetscFunctionReturn(0);
637b5b4017aSHong Zhang }
638b5b4017aSHong Zhang 
639b5b4017aSHong Zhang /*@C
640*b1e111ebSHong Zhang   TSComputeIHessianProductFunctionPP - Runs the user-defined vector-Hessian-vector product function for Fpp.
641b5b4017aSHong Zhang 
642b5b4017aSHong Zhang   Collective on TS
643b5b4017aSHong Zhang 
644b5b4017aSHong Zhang   Input Parameters:
645b5b4017aSHong Zhang . ts   - The TS context obtained from TSCreate()
646b5b4017aSHong Zhang 
647b5b4017aSHong Zhang   Notes:
648*b1e111ebSHong Zhang   TSComputeIHessianProductFunctionPP() is typically used for sensitivity implementation,
649b5b4017aSHong Zhang   so most users would not generally call this routine themselves.
650b5b4017aSHong Zhang 
651b5b4017aSHong Zhang   Level: developer
652b5b4017aSHong Zhang 
653b5b4017aSHong Zhang .keywords: TS, sensitivity
654b5b4017aSHong Zhang 
655b5b4017aSHong Zhang .seealso: TSSetIHessianProduct()
656b5b4017aSHong Zhang @*/
657*b1e111ebSHong Zhang PetscErrorCode TSComputeIHessianProductFunctionPP(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV)
658b5b4017aSHong Zhang {
659b5b4017aSHong Zhang   PetscErrorCode ierr;
660b5b4017aSHong Zhang 
661b5b4017aSHong Zhang   PetscFunctionBegin;
66233f72c5dSHong Zhang   if (!VHV) PetscFunctionReturn(0);
663b5b4017aSHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
664b5b4017aSHong Zhang   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
665b5b4017aSHong Zhang 
66667633408SHong Zhang   if (ts->ihessianproduct_fpp) {
667b5b4017aSHong Zhang     PetscStackPush("TS user IHessianProduct function 3 for sensitivity analysis");
668b5b4017aSHong Zhang     ierr = (*ts->ihessianproduct_fpp)(ts,t,U,Vl,Vr,VHV,ts->ihessianproductctx);CHKERRQ(ierr);
669b5b4017aSHong Zhang     PetscStackPop;
67067633408SHong Zhang   }
67167633408SHong Zhang   /* does not consider IMEX for now, so either IHessian or RHSHessian will be calculated, using the same output VHV */
67267633408SHong Zhang   if (ts->rhshessianproduct_gpp) {
67367633408SHong Zhang     PetscInt nadj;
674*b1e111ebSHong Zhang     ierr = TSComputeRHSHessianProductFunctionPP(ts,t,U,Vl,Vr,VHV);CHKERRQ(ierr);
67567633408SHong Zhang     for (nadj=0; nadj<ts->numcost; nadj++) {
67667633408SHong Zhang       ierr = VecScale(VHV[nadj],-1);CHKERRQ(ierr);
67767633408SHong Zhang     }
67867633408SHong Zhang   }
679b5b4017aSHong Zhang   PetscFunctionReturn(0);
680b5b4017aSHong Zhang }
681b5b4017aSHong Zhang 
68213af1a74SHong Zhang /*@C
68313af1a74SHong Zhang   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.
68413af1a74SHong Zhang 
68513af1a74SHong Zhang   Logically Collective on TS
68613af1a74SHong Zhang 
68713af1a74SHong Zhang   Input Parameters:
68813af1a74SHong Zhang + ts - TS context obtained from TSCreate()
68913af1a74SHong Zhang . rhshp1 - an array of vectors storing the result of vector-Hessian-vector product for G_UU
69013af1a74SHong Zhang . hessianproductfunc1 - vector-Hessian-vector product function for G_UU
69113af1a74SHong Zhang . rhshp2 - an array of vectors storing the result of vector-Hessian-vector product for G_UP
69213af1a74SHong Zhang . hessianproductfunc2 - vector-Hessian-vector product function for G_UP
69313af1a74SHong Zhang . rhshp3 - an array of vectors storing the result of vector-Hessian-vector product for G_PU
69413af1a74SHong Zhang . hessianproductfunc3 - vector-Hessian-vector product function for G_PU
69513af1a74SHong Zhang . rhshp4 - an array of vectors storing the result of vector-Hessian-vector product for G_PP
69613af1a74SHong Zhang . hessianproductfunc4 - vector-Hessian-vector product function for G_PP
69713af1a74SHong Zhang 
69813af1a74SHong Zhang   Calling sequence of ihessianproductfunc:
699369cf35fSHong Zhang $ rhshessianproductfunc (TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx);
70013af1a74SHong Zhang +   t - current timestep
70113af1a74SHong Zhang .   U - input vector (current ODE solution)
702369cf35fSHong Zhang .   Vl - an array of input vectors to be left-multiplied with the Hessian
70313af1a74SHong Zhang .   Vr - input vector to be right-multiplied with the Hessian
704369cf35fSHong Zhang .   VHV - an array of output vectors for vector-Hessian-vector product
70513af1a74SHong Zhang -   ctx - [optional] user-defined function context
70613af1a74SHong Zhang 
70713af1a74SHong Zhang   Level: intermediate
70813af1a74SHong Zhang 
709369cf35fSHong Zhang   Notes:
710369cf35fSHong Zhang   The first Hessian function and the working array are required.
711369cf35fSHong Zhang   As an example to implement the callback functions, the second callback function calculates the vector-Hessian-vector product
712369cf35fSHong Zhang   $ Vl_n^T*G_UP*Vr
713369cf35fSHong Zhang   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.
714369cf35fSHong Zhang   Each entry of G_UP corresponds to the derivative
715369cf35fSHong Zhang   $ G_UP[i][j][k] = \frac{\partial^2 G[i]}{\partial U[j] \partial P[k]}.
716369cf35fSHong Zhang   The result of the vector-Hessian-vector product for Vl_n needs to be stored in vector VHV_n with j-th entry being
717369cf35fSHong Zhang   $ VHV_n[j] = \sum_i \sum_k {Vl_n[i] * G_UP[i][j][k] * Vr[k]}
718369cf35fSHong Zhang   If the cost function is a scalar, there will be only one vector in Vl and VHV.
71913af1a74SHong Zhang 
72013af1a74SHong Zhang .keywords: TS, sensitivity
72113af1a74SHong Zhang 
72213af1a74SHong Zhang .seealso:
72313af1a74SHong Zhang @*/
72413af1a74SHong Zhang PetscErrorCode TSSetRHSHessianProduct(TS ts,Vec *rhshp1,PetscErrorCode (*rhshessianproductfunc1)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*),
72513af1a74SHong Zhang                                           Vec *rhshp2,PetscErrorCode (*rhshessianproductfunc2)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*),
72613af1a74SHong Zhang                                           Vec *rhshp3,PetscErrorCode (*rhshessianproductfunc3)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*),
72713af1a74SHong Zhang                                           Vec *rhshp4,PetscErrorCode (*rhshessianproductfunc4)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*),
72813af1a74SHong Zhang                                     void *ctx)
72913af1a74SHong Zhang {
73013af1a74SHong Zhang   PetscFunctionBegin;
73113af1a74SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
73213af1a74SHong Zhang   PetscValidPointer(rhshp1,2);
73313af1a74SHong Zhang 
73413af1a74SHong Zhang   ts->rhshessianproductctx = ctx;
73513af1a74SHong Zhang   if (rhshp1) ts->vecs_guu = rhshp1;
73613af1a74SHong Zhang   if (rhshp2) ts->vecs_gup = rhshp2;
73713af1a74SHong Zhang   if (rhshp3) ts->vecs_gpu = rhshp3;
73813af1a74SHong Zhang   if (rhshp4) ts->vecs_gpp = rhshp4;
73913af1a74SHong Zhang   ts->rhshessianproduct_guu = rhshessianproductfunc1;
74013af1a74SHong Zhang   ts->rhshessianproduct_gup = rhshessianproductfunc2;
74113af1a74SHong Zhang   ts->rhshessianproduct_gpu = rhshessianproductfunc3;
74213af1a74SHong Zhang   ts->rhshessianproduct_gpp = rhshessianproductfunc4;
74313af1a74SHong Zhang   PetscFunctionReturn(0);
74413af1a74SHong Zhang }
74513af1a74SHong Zhang 
74613af1a74SHong Zhang /*@C
747*b1e111ebSHong Zhang   TSComputeRHSHessianProductFunctionUU - Runs the user-defined vector-Hessian-vector product function for Guu.
74813af1a74SHong Zhang 
74913af1a74SHong Zhang   Collective on TS
75013af1a74SHong Zhang 
75113af1a74SHong Zhang   Input Parameters:
75213af1a74SHong Zhang . ts   - The TS context obtained from TSCreate()
75313af1a74SHong Zhang 
75413af1a74SHong Zhang   Notes:
755*b1e111ebSHong Zhang   TSComputeRHSHessianProductFunctionUU() is typically used for sensitivity implementation,
75613af1a74SHong Zhang   so most users would not generally call this routine themselves.
75713af1a74SHong Zhang 
75813af1a74SHong Zhang   Level: developer
75913af1a74SHong Zhang 
76013af1a74SHong Zhang .keywords: TS, sensitivity
76113af1a74SHong Zhang 
76213af1a74SHong Zhang .seealso: TSSetRHSHessianProduct()
76313af1a74SHong Zhang @*/
764*b1e111ebSHong Zhang PetscErrorCode TSComputeRHSHessianProductFunctionUU(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV)
76513af1a74SHong Zhang {
76613af1a74SHong Zhang   PetscErrorCode ierr;
76713af1a74SHong Zhang 
76813af1a74SHong Zhang   PetscFunctionBegin;
76913af1a74SHong Zhang   if (!VHV) PetscFunctionReturn(0);
77013af1a74SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
77113af1a74SHong Zhang   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
77213af1a74SHong Zhang 
77313af1a74SHong Zhang   PetscStackPush("TS user RHSHessianProduct function 1 for sensitivity analysis");
77413af1a74SHong Zhang   ierr = (*ts->rhshessianproduct_guu)(ts,t,U,Vl,Vr,VHV,ts->rhshessianproductctx);CHKERRQ(ierr);
77513af1a74SHong Zhang   PetscStackPop;
77613af1a74SHong Zhang   PetscFunctionReturn(0);
77713af1a74SHong Zhang }
77813af1a74SHong Zhang 
77913af1a74SHong Zhang /*@C
780*b1e111ebSHong Zhang   TSComputeRHSHessianProductFunctionUP - Runs the user-defined vector-Hessian-vector product function for Gup.
78113af1a74SHong Zhang 
78213af1a74SHong Zhang   Collective on TS
78313af1a74SHong Zhang 
78413af1a74SHong Zhang   Input Parameters:
78513af1a74SHong Zhang . ts   - The TS context obtained from TSCreate()
78613af1a74SHong Zhang 
78713af1a74SHong Zhang   Notes:
788*b1e111ebSHong Zhang   TSComputeRHSHessianProductFunctionUP() is typically used for sensitivity implementation,
78913af1a74SHong Zhang   so most users would not generally call this routine themselves.
79013af1a74SHong Zhang 
79113af1a74SHong Zhang   Level: developer
79213af1a74SHong Zhang 
79313af1a74SHong Zhang .keywords: TS, sensitivity
79413af1a74SHong Zhang 
79513af1a74SHong Zhang .seealso: TSSetRHSHessianProduct()
79613af1a74SHong Zhang @*/
797*b1e111ebSHong Zhang PetscErrorCode TSComputeRHSHessianProductFunctionUP(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV)
79813af1a74SHong Zhang {
79913af1a74SHong Zhang   PetscErrorCode ierr;
80013af1a74SHong Zhang 
80113af1a74SHong Zhang   PetscFunctionBegin;
80213af1a74SHong Zhang   if (!VHV) PetscFunctionReturn(0);
80313af1a74SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
80413af1a74SHong Zhang   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
80513af1a74SHong Zhang 
80613af1a74SHong Zhang   PetscStackPush("TS user RHSHessianProduct function 2 for sensitivity analysis");
80713af1a74SHong Zhang   ierr = (*ts->rhshessianproduct_gup)(ts,t,U,Vl,Vr,VHV,ts->rhshessianproductctx);CHKERRQ(ierr);
80813af1a74SHong Zhang   PetscStackPop;
80913af1a74SHong Zhang   PetscFunctionReturn(0);
81013af1a74SHong Zhang }
81113af1a74SHong Zhang 
81213af1a74SHong Zhang /*@C
813*b1e111ebSHong Zhang   TSComputeRHSHessianProductFunctionPU - Runs the user-defined vector-Hessian-vector product function for Gpu.
81413af1a74SHong Zhang 
81513af1a74SHong Zhang   Collective on TS
81613af1a74SHong Zhang 
81713af1a74SHong Zhang   Input Parameters:
81813af1a74SHong Zhang . ts   - The TS context obtained from TSCreate()
81913af1a74SHong Zhang 
82013af1a74SHong Zhang   Notes:
821*b1e111ebSHong Zhang   TSComputeRHSHessianProductFunctionPU() is typically used for sensitivity implementation,
82213af1a74SHong Zhang   so most users would not generally call this routine themselves.
82313af1a74SHong Zhang 
82413af1a74SHong Zhang   Level: developer
82513af1a74SHong Zhang 
82613af1a74SHong Zhang .keywords: TS, sensitivity
82713af1a74SHong Zhang 
82813af1a74SHong Zhang .seealso: TSSetRHSHessianProduct()
82913af1a74SHong Zhang @*/
830*b1e111ebSHong Zhang PetscErrorCode TSComputeRHSHessianProductFunctionPU(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV)
83113af1a74SHong Zhang {
83213af1a74SHong Zhang   PetscErrorCode ierr;
83313af1a74SHong Zhang 
83413af1a74SHong Zhang   PetscFunctionBegin;
83513af1a74SHong Zhang   if (!VHV) PetscFunctionReturn(0);
83613af1a74SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
83713af1a74SHong Zhang   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
83813af1a74SHong Zhang 
83913af1a74SHong Zhang   PetscStackPush("TS user RHSHessianProduct function 3 for sensitivity analysis");
84013af1a74SHong Zhang   ierr = (*ts->rhshessianproduct_gpu)(ts,t,U,Vl,Vr,VHV,ts->rhshessianproductctx);CHKERRQ(ierr);
84113af1a74SHong Zhang   PetscStackPop;
84213af1a74SHong Zhang   PetscFunctionReturn(0);
84313af1a74SHong Zhang }
84413af1a74SHong Zhang 
84513af1a74SHong Zhang /*@C
846*b1e111ebSHong Zhang   TSComputeRHSHessianProductFunctionPP - Runs the user-defined vector-Hessian-vector product function for Gpp.
84713af1a74SHong Zhang 
84813af1a74SHong Zhang   Collective on TS
84913af1a74SHong Zhang 
85013af1a74SHong Zhang   Input Parameters:
85113af1a74SHong Zhang . ts   - The TS context obtained from TSCreate()
85213af1a74SHong Zhang 
85313af1a74SHong Zhang   Notes:
854*b1e111ebSHong Zhang   TSComputeRHSHessianProductFunctionPP() is typically used for sensitivity implementation,
85513af1a74SHong Zhang   so most users would not generally call this routine themselves.
85613af1a74SHong Zhang 
85713af1a74SHong Zhang   Level: developer
85813af1a74SHong Zhang 
85913af1a74SHong Zhang .keywords: TS, sensitivity
86013af1a74SHong Zhang 
86113af1a74SHong Zhang .seealso: TSSetRHSHessianProduct()
86213af1a74SHong Zhang @*/
863*b1e111ebSHong Zhang PetscErrorCode TSComputeRHSHessianProductFunctionPP(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV)
86413af1a74SHong Zhang {
86513af1a74SHong Zhang   PetscErrorCode ierr;
86613af1a74SHong Zhang 
86713af1a74SHong Zhang   PetscFunctionBegin;
86813af1a74SHong Zhang   if (!VHV) PetscFunctionReturn(0);
86913af1a74SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
87013af1a74SHong Zhang   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
87113af1a74SHong Zhang 
87213af1a74SHong Zhang   PetscStackPush("TS user RHSHessianProduct function 3 for sensitivity analysis");
87313af1a74SHong Zhang   ierr = (*ts->rhshessianproduct_gpp)(ts,t,U,Vl,Vr,VHV,ts->rhshessianproductctx);CHKERRQ(ierr);
87413af1a74SHong Zhang   PetscStackPop;
87513af1a74SHong Zhang   PetscFunctionReturn(0);
87613af1a74SHong Zhang }
87713af1a74SHong Zhang 
878814e85d6SHong Zhang /* --------------------------- Adjoint sensitivity ---------------------------*/
879814e85d6SHong Zhang 
880814e85d6SHong Zhang /*@
881814e85d6SHong Zhang    TSSetCostGradients - Sets the initial value of the gradients of the cost function w.r.t. initial values and w.r.t. the problem parameters
882814e85d6SHong Zhang       for use by the TSAdjoint routines.
883814e85d6SHong Zhang 
884814e85d6SHong Zhang    Logically Collective on TS and Vec
885814e85d6SHong Zhang 
886814e85d6SHong Zhang    Input Parameters:
887814e85d6SHong Zhang +  ts - the TS context obtained from TSCreate()
888814e85d6SHong Zhang .  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
889814e85d6SHong Zhang -  mu - gradients with respect to the parameters, the number of entries in these vectors is the same as the number of parameters
890814e85d6SHong Zhang 
891814e85d6SHong Zhang    Level: beginner
892814e85d6SHong Zhang 
893814e85d6SHong Zhang    Notes:
894814e85d6SHong Zhang     the entries in these vectors must be correctly initialized with the values lamda_i = df/dy|finaltime  mu_i = df/dp|finaltime
895814e85d6SHong Zhang 
896814e85d6SHong Zhang    After TSAdjointSolve() is called the lamba and the mu contain the computed sensitivities
897814e85d6SHong Zhang 
898814e85d6SHong Zhang .keywords: TS, timestep, set, sensitivity, initial values
899b5b4017aSHong Zhang 
900b5b4017aSHong Zhang .seealso TSGetCostGradients()
901814e85d6SHong Zhang @*/
902814e85d6SHong Zhang PetscErrorCode TSSetCostGradients(TS ts,PetscInt numcost,Vec *lambda,Vec *mu)
903814e85d6SHong Zhang {
904814e85d6SHong Zhang   PetscFunctionBegin;
905814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
906814e85d6SHong Zhang   PetscValidPointer(lambda,2);
907814e85d6SHong Zhang   ts->vecs_sensi  = lambda;
908814e85d6SHong Zhang   ts->vecs_sensip = mu;
909814e85d6SHong Zhang   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");
910814e85d6SHong Zhang   ts->numcost  = numcost;
911814e85d6SHong Zhang   PetscFunctionReturn(0);
912814e85d6SHong Zhang }
913814e85d6SHong Zhang 
914814e85d6SHong Zhang /*@
915814e85d6SHong Zhang    TSGetCostGradients - Returns the gradients from the TSAdjointSolve()
916814e85d6SHong Zhang 
917814e85d6SHong Zhang    Not Collective, but Vec returned is parallel if TS is parallel
918814e85d6SHong Zhang 
919814e85d6SHong Zhang    Input Parameter:
920814e85d6SHong Zhang .  ts - the TS context obtained from TSCreate()
921814e85d6SHong Zhang 
922814e85d6SHong Zhang    Output Parameter:
923814e85d6SHong Zhang +  lambda - vectors containing the gradients of the cost functions with respect to the ODE/DAE solution variables
924814e85d6SHong Zhang -  mu - vectors containing the gradients of the cost functions with respect to the problem parameters
925814e85d6SHong Zhang 
926814e85d6SHong Zhang    Level: intermediate
927814e85d6SHong Zhang 
928814e85d6SHong Zhang .keywords: TS, timestep, get, sensitivity
929b5b4017aSHong Zhang 
930b5b4017aSHong Zhang .seealso: TSSetCostGradients()
931814e85d6SHong Zhang @*/
932814e85d6SHong Zhang PetscErrorCode TSGetCostGradients(TS ts,PetscInt *numcost,Vec **lambda,Vec **mu)
933814e85d6SHong Zhang {
934814e85d6SHong Zhang   PetscFunctionBegin;
935814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
936814e85d6SHong Zhang   if (numcost) *numcost = ts->numcost;
937814e85d6SHong Zhang   if (lambda)  *lambda  = ts->vecs_sensi;
938814e85d6SHong Zhang   if (mu)      *mu      = ts->vecs_sensip;
939814e85d6SHong Zhang   PetscFunctionReturn(0);
940814e85d6SHong Zhang }
941814e85d6SHong Zhang 
942814e85d6SHong Zhang /*@
943b5b4017aSHong Zhang    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
944b5b4017aSHong Zhang       for use by the TSAdjoint routines.
945b5b4017aSHong Zhang 
946b5b4017aSHong Zhang    Logically Collective on TS and Vec
947b5b4017aSHong Zhang 
948b5b4017aSHong Zhang    Input Parameters:
949b5b4017aSHong Zhang +  ts - the TS context obtained from TSCreate()
950b5b4017aSHong Zhang .  numcost - number of cost functions
951b5b4017aSHong Zhang .  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
952b5b4017aSHong Zhang .  mu2 - Hessian-vector product with respect to the parameters, the number of entries in these vectors is the same as the number of parameters
953b5b4017aSHong Zhang -  dir - the direction vector that are multiplied with the Hessian of the cost functions
954b5b4017aSHong Zhang 
955b5b4017aSHong Zhang    Level: beginner
956b5b4017aSHong Zhang 
957b5b4017aSHong Zhang    Notes: Hessian of the cost function is completely different from Hessian of the ODE/DAE system
958b5b4017aSHong Zhang 
959ecf68647SHong Zhang    For second-order adjoint, one needs to call this function and then TSAdjointSetForward() before TSSolve().
960b5b4017aSHong Zhang 
961b5b4017aSHong Zhang    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.
962b5b4017aSHong Zhang 
9633fca9621SHong Zhang    Passing NULL for lambda2 disables the second-order calculation.
964b5b4017aSHong Zhang .keywords: TS, sensitivity, second-order adjoint
965b5b4017aSHong Zhang 
966ecf68647SHong Zhang .seealso: TSAdjointSetForward()
967b5b4017aSHong Zhang @*/
968b5b4017aSHong Zhang PetscErrorCode TSSetCostHessianProducts(TS ts,PetscInt numcost,Vec *lambda2,Vec *mu2,Vec dir)
969b5b4017aSHong Zhang {
970b5b4017aSHong Zhang   PetscFunctionBegin;
971b5b4017aSHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
972b5b4017aSHong Zhang   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");
973b5b4017aSHong Zhang   ts->numcost       = numcost;
974b5b4017aSHong Zhang   ts->vecs_sensi2   = lambda2;
97533f72c5dSHong Zhang   ts->vecs_sensi2p  = mu2;
976b5b4017aSHong Zhang   ts->vec_dir       = dir;
977b5b4017aSHong Zhang   PetscFunctionReturn(0);
978b5b4017aSHong Zhang }
979b5b4017aSHong Zhang 
980b5b4017aSHong Zhang /*@
981b5b4017aSHong Zhang    TSGetCostHessianProducts - Returns the gradients from the TSAdjointSolve()
982b5b4017aSHong Zhang 
983b5b4017aSHong Zhang    Not Collective, but Vec returned is parallel if TS is parallel
984b5b4017aSHong Zhang 
985b5b4017aSHong Zhang    Input Parameter:
986b5b4017aSHong Zhang .  ts - the TS context obtained from TSCreate()
987b5b4017aSHong Zhang 
988b5b4017aSHong Zhang    Output Parameter:
989b5b4017aSHong Zhang +  numcost - number of cost functions
990b5b4017aSHong Zhang .  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
991b5b4017aSHong Zhang .  mu2 - Hessian-vector product with respect to the parameters, the number of entries in these vectors is the same as the number of parameters
992b5b4017aSHong Zhang -  dir - the direction vector that are multiplied with the Hessian of the cost functions
993b5b4017aSHong Zhang 
994b5b4017aSHong Zhang    Level: intermediate
995b5b4017aSHong Zhang 
996b5b4017aSHong Zhang .keywords: TS, sensitivity, second-order adjoint
997b5b4017aSHong Zhang 
998b5b4017aSHong Zhang .seealso: TSSetCostHessianProducts()
999b5b4017aSHong Zhang @*/
1000b5b4017aSHong Zhang PetscErrorCode TSGetCostHessianProducts(TS ts,PetscInt *numcost,Vec **lambda2,Vec **mu2, Vec *dir)
1001b5b4017aSHong Zhang {
1002b5b4017aSHong Zhang   PetscFunctionBegin;
1003b5b4017aSHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1004b5b4017aSHong Zhang   if (numcost) *numcost = ts->numcost;
1005b5b4017aSHong Zhang   if (lambda2) *lambda2 = ts->vecs_sensi2;
100633f72c5dSHong Zhang   if (mu2)     *mu2     = ts->vecs_sensi2p;
1007b5b4017aSHong Zhang   if (dir)     *dir     = ts->vec_dir;
1008b5b4017aSHong Zhang   PetscFunctionReturn(0);
1009b5b4017aSHong Zhang }
1010b5b4017aSHong Zhang 
1011b5b4017aSHong Zhang /*@
1012ecf68647SHong Zhang   TSAdjointSetForward - Trigger the tangent linear solver and initialize the forward sensitivities
1013b5b4017aSHong Zhang 
1014b5b4017aSHong Zhang   Logically Collective on TS and Mat
1015b5b4017aSHong Zhang 
1016b5b4017aSHong Zhang   Input Parameters:
1017b5b4017aSHong Zhang +  ts - the TS context obtained from TSCreate()
1018b5b4017aSHong Zhang -  didp - the derivative of initial values w.r.t. parameters
1019b5b4017aSHong Zhang 
1020b5b4017aSHong Zhang   Level: intermediate
1021b5b4017aSHong Zhang 
1022ecf68647SHong Zhang   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.
1023b5b4017aSHong Zhang 
1024b5b4017aSHong Zhang .keywords: TS, sensitivity, second-order adjoint
1025b5b4017aSHong Zhang 
1026ecf68647SHong Zhang .seealso: TSSetCostHessianProducts(), TSAdjointResetForward()
1027b5b4017aSHong Zhang @*/
1028ecf68647SHong Zhang PetscErrorCode TSAdjointSetForward(TS ts,Mat didp)
1029b5b4017aSHong Zhang {
103033f72c5dSHong Zhang   Mat            A;
103133f72c5dSHong Zhang   Vec            sp;
103233f72c5dSHong Zhang   PetscScalar    *xarr;
103333f72c5dSHong Zhang   PetscInt       lsize;
1034b5b4017aSHong Zhang   PetscErrorCode ierr;
1035b5b4017aSHong Zhang 
1036b5b4017aSHong Zhang   PetscFunctionBegin;
1037b5b4017aSHong Zhang   ts->forward_solve = PETSC_TRUE; /* turn on tangent linear mode */
1038b5b4017aSHong Zhang   if (!ts->vecs_sensi2) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must call TSSetCostHessianProducts() first");
103933f72c5dSHong Zhang   if (!ts->vec_dir) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Directional vector is missing. Call TSSetCostHessianProducts() to set it.");
104033f72c5dSHong Zhang   /* create a single-column dense matrix */
104133f72c5dSHong Zhang   ierr = VecGetLocalSize(ts->vec_sol,&lsize);CHKERRQ(ierr);
1042f63bf25fSHong Zhang   ierr = MatCreateDense(PetscObjectComm((PetscObject)ts),lsize,PETSC_DECIDE,PETSC_DECIDE,1,NULL,&A);CHKERRQ(ierr);
104333f72c5dSHong Zhang 
104433f72c5dSHong Zhang   ierr = VecDuplicate(ts->vec_sol,&sp);CHKERRQ(ierr);
104533f72c5dSHong Zhang   ierr = MatDenseGetColumn(A,0,&xarr);CHKERRQ(ierr);
104633f72c5dSHong Zhang   ierr = VecPlaceArray(sp,xarr);CHKERRQ(ierr);
1047ecf68647SHong Zhang   if (ts->vecs_sensi2p) { /* tangent linear variable initialized as 2*dIdP*dir */
104833f72c5dSHong Zhang     if (didp) {
104933f72c5dSHong Zhang       ierr = MatMult(didp,ts->vec_dir,sp);CHKERRQ(ierr);
105033f72c5dSHong Zhang       ierr = VecScale(sp,2.);CHKERRQ(ierr);
105133f72c5dSHong Zhang     } else {
105233f72c5dSHong Zhang       ierr = VecZeroEntries(sp);CHKERRQ(ierr);
105333f72c5dSHong Zhang     }
1054ecf68647SHong Zhang   } else { /* tangent linear variable initialized as dir */
105533f72c5dSHong Zhang     ierr = VecCopy(ts->vec_dir,sp);CHKERRQ(ierr);
105633f72c5dSHong Zhang   }
105733f72c5dSHong Zhang   ierr = VecResetArray(sp);CHKERRQ(ierr);
105833f72c5dSHong Zhang   ierr = MatDenseRestoreColumn(A,&xarr);CHKERRQ(ierr);
105933f72c5dSHong Zhang   ierr = VecDestroy(&sp);CHKERRQ(ierr);
106033f72c5dSHong Zhang 
106133f72c5dSHong Zhang   ierr = TSForwardSetInitialSensitivities(ts,A);CHKERRQ(ierr); /* if didp is NULL, identity matrix is assumed */
106233f72c5dSHong Zhang 
106333f72c5dSHong Zhang   ierr = MatDestroy(&A);CHKERRQ(ierr);
1064b5b4017aSHong Zhang   PetscFunctionReturn(0);
1065b5b4017aSHong Zhang }
1066b5b4017aSHong Zhang 
1067b5b4017aSHong Zhang /*@
1068ecf68647SHong Zhang   TSAdjointResetForward - Reset the tangent linear solver and destroy the tangent linear context
1069ecf68647SHong Zhang 
1070ecf68647SHong Zhang   Logically Collective on TS and Mat
1071ecf68647SHong Zhang 
1072ecf68647SHong Zhang   Input Parameters:
1073ecf68647SHong Zhang .  ts - the TS context obtained from TSCreate()
1074ecf68647SHong Zhang 
1075ecf68647SHong Zhang   Level: intermediate
1076ecf68647SHong Zhang 
1077ecf68647SHong Zhang .keywords: TS, sensitivity, second-order adjoint
1078ecf68647SHong Zhang 
1079ecf68647SHong Zhang .seealso: TSAdjointSetForward()
1080ecf68647SHong Zhang @*/
1081ecf68647SHong Zhang PetscErrorCode TSAdjointResetForward(TS ts)
1082ecf68647SHong Zhang {
1083ecf68647SHong Zhang   PetscErrorCode ierr;
1084ecf68647SHong Zhang 
1085ecf68647SHong Zhang   PetscFunctionBegin;
1086ecf68647SHong Zhang   ts->forward_solve = PETSC_FALSE; /* turn off tangent linear mode */
1087ecf68647SHong Zhang   ierr = TSForwardReset(ts);CHKERRQ(ierr);
1088ecf68647SHong Zhang   PetscFunctionReturn(0);
1089ecf68647SHong Zhang }
1090ecf68647SHong Zhang 
1091ecf68647SHong Zhang /*@
1092814e85d6SHong Zhang    TSAdjointSetUp - Sets up the internal data structures for the later use
1093814e85d6SHong Zhang    of an adjoint solver
1094814e85d6SHong Zhang 
1095814e85d6SHong Zhang    Collective on TS
1096814e85d6SHong Zhang 
1097814e85d6SHong Zhang    Input Parameter:
1098814e85d6SHong Zhang .  ts - the TS context obtained from TSCreate()
1099814e85d6SHong Zhang 
1100814e85d6SHong Zhang    Level: advanced
1101814e85d6SHong Zhang 
1102814e85d6SHong Zhang .keywords: TS, timestep, setup
1103814e85d6SHong Zhang 
1104814e85d6SHong Zhang .seealso: TSCreate(), TSAdjointStep(), TSSetCostGradients()
1105814e85d6SHong Zhang @*/
1106814e85d6SHong Zhang PetscErrorCode TSAdjointSetUp(TS ts)
1107814e85d6SHong Zhang {
1108881c1a9bSHong Zhang   TSTrajectory     tj;
1109881c1a9bSHong Zhang   PetscBool        match;
1110814e85d6SHong Zhang   PetscErrorCode   ierr;
1111814e85d6SHong Zhang 
1112814e85d6SHong Zhang   PetscFunctionBegin;
1113814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1114814e85d6SHong Zhang   if (ts->adjointsetupcalled) PetscFunctionReturn(0);
1115814e85d6SHong Zhang   if (!ts->vecs_sensi) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetCostGradients() first");
111633f72c5dSHong Zhang   if (ts->vecs_sensip && !ts->Jacp && !ts->Jacprhs) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetRHSJacobianP() or TSSetIJacobianP() first");
1117881c1a9bSHong Zhang   ierr = TSGetTrajectory(ts,&tj);CHKERRQ(ierr);
11188e224c67SHong Zhang   ierr = PetscObjectTypeCompare((PetscObject)tj,TSTRAJECTORYBASIC,&match);CHKERRQ(ierr);
1119881c1a9bSHong Zhang   if (match) {
1120881c1a9bSHong Zhang     PetscBool solution_only;
1121881c1a9bSHong Zhang     ierr = TSTrajectoryGetSolutionOnly(tj,&solution_only);CHKERRQ(ierr);
1122881c1a9bSHong Zhang     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");
1123881c1a9bSHong Zhang   }
112433f72c5dSHong Zhang 
112533f72c5dSHong Zhang   if (!ts->Jacp && ts->Jacprhs) ts->Jacp = ts->Jacprhs;
1126814e85d6SHong Zhang 
1127cd4cee2dSHong Zhang   if (ts->quadraturets) { /* if there is integral in the cost function */
1128cd4cee2dSHong Zhang     ierr = VecDuplicate(ts->vecs_sensi[0],&ts->vec_drdu_col);CHKERRQ(ierr);
1129814e85d6SHong Zhang     if (ts->vecs_sensip){
1130cd4cee2dSHong Zhang       ierr = VecDuplicate(ts->vecs_sensip[0],&ts->vec_drdp_col);CHKERRQ(ierr);
1131814e85d6SHong Zhang     }
1132814e85d6SHong Zhang   }
1133814e85d6SHong Zhang 
1134814e85d6SHong Zhang   if (ts->ops->adjointsetup) {
1135814e85d6SHong Zhang     ierr = (*ts->ops->adjointsetup)(ts);CHKERRQ(ierr);
1136814e85d6SHong Zhang   }
1137814e85d6SHong Zhang   ts->adjointsetupcalled = PETSC_TRUE;
1138814e85d6SHong Zhang   PetscFunctionReturn(0);
1139814e85d6SHong Zhang }
1140814e85d6SHong Zhang 
1141814e85d6SHong Zhang /*@
1142ece46509SHong Zhang    TSAdjointReset - Resets a TSAdjoint context and removes any allocated Vecs and Mats.
1143ece46509SHong Zhang 
1144ece46509SHong Zhang    Collective on TS
1145ece46509SHong Zhang 
1146ece46509SHong Zhang    Input Parameter:
1147ece46509SHong Zhang .  ts - the TS context obtained from TSCreate()
1148ece46509SHong Zhang 
1149ece46509SHong Zhang    Level: beginner
1150ece46509SHong Zhang 
1151ece46509SHong Zhang .keywords: TS, timestep, reset
1152ece46509SHong Zhang 
1153ece46509SHong Zhang .seealso: TSCreate(), TSAdjointSetup(), TSADestroy()
1154ece46509SHong Zhang @*/
1155ece46509SHong Zhang PetscErrorCode TSAdjointReset(TS ts)
1156ece46509SHong Zhang {
1157ece46509SHong Zhang   PetscErrorCode ierr;
1158ece46509SHong Zhang 
1159ece46509SHong Zhang   PetscFunctionBegin;
1160ece46509SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1161ece46509SHong Zhang   if (ts->ops->adjointreset) {
1162ece46509SHong Zhang     ierr = (*ts->ops->adjointreset)(ts);CHKERRQ(ierr);
1163ece46509SHong Zhang   }
11647207e0fdSHong Zhang   if (ts->quadraturets) { /* if there is integral in the cost function */
11657207e0fdSHong Zhang     ierr = VecDestroy(&ts->vec_drdu_col);CHKERRQ(ierr);
11667207e0fdSHong Zhang     if (ts->vecs_sensip){
11677207e0fdSHong Zhang       ierr = VecDestroy(&ts->vec_drdp_col);CHKERRQ(ierr);
11687207e0fdSHong Zhang     }
11697207e0fdSHong Zhang   }
1170ece46509SHong Zhang   ts->vecs_sensi         = NULL;
1171ece46509SHong Zhang   ts->vecs_sensip        = NULL;
1172ece46509SHong Zhang   ts->vecs_sensi2        = NULL;
117333f72c5dSHong Zhang   ts->vecs_sensi2p       = NULL;
1174ece46509SHong Zhang   ts->vec_dir            = NULL;
1175ece46509SHong Zhang   ts->adjointsetupcalled = PETSC_FALSE;
1176ece46509SHong Zhang   PetscFunctionReturn(0);
1177ece46509SHong Zhang }
1178ece46509SHong Zhang 
1179ece46509SHong Zhang /*@
1180814e85d6SHong Zhang    TSAdjointSetSteps - Sets the number of steps the adjoint solver should take backward in time
1181814e85d6SHong Zhang 
1182814e85d6SHong Zhang    Logically Collective on TS
1183814e85d6SHong Zhang 
1184814e85d6SHong Zhang    Input Parameters:
1185814e85d6SHong Zhang +  ts - the TS context obtained from TSCreate()
1186814e85d6SHong Zhang .  steps - number of steps to use
1187814e85d6SHong Zhang 
1188814e85d6SHong Zhang    Level: intermediate
1189814e85d6SHong Zhang 
1190814e85d6SHong Zhang    Notes:
1191814e85d6SHong Zhang     Normally one does not call this and TSAdjointSolve() integrates back to the original timestep. One can call this
1192814e85d6SHong Zhang           so as to integrate back to less than the original timestep
1193814e85d6SHong Zhang 
1194814e85d6SHong Zhang .keywords: TS, timestep, set, maximum, iterations
1195814e85d6SHong Zhang 
1196814e85d6SHong Zhang .seealso: TSSetExactFinalTime()
1197814e85d6SHong Zhang @*/
1198814e85d6SHong Zhang PetscErrorCode TSAdjointSetSteps(TS ts,PetscInt steps)
1199814e85d6SHong Zhang {
1200814e85d6SHong Zhang   PetscFunctionBegin;
1201814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1202814e85d6SHong Zhang   PetscValidLogicalCollectiveInt(ts,steps,2);
1203814e85d6SHong Zhang   if (steps < 0) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Cannot step back a negative number of steps");
1204814e85d6SHong Zhang   if (steps > ts->steps) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Cannot step back more than the total number of forward steps");
1205814e85d6SHong Zhang   ts->adjoint_max_steps = steps;
1206814e85d6SHong Zhang   PetscFunctionReturn(0);
1207814e85d6SHong Zhang }
1208814e85d6SHong Zhang 
1209814e85d6SHong Zhang /*@C
1210814e85d6SHong Zhang   TSAdjointSetRHSJacobian - Deprecated, use TSSetRHSJacobianP()
1211814e85d6SHong Zhang 
1212814e85d6SHong Zhang   Level: deprecated
1213814e85d6SHong Zhang 
1214814e85d6SHong Zhang @*/
1215814e85d6SHong Zhang PetscErrorCode TSAdjointSetRHSJacobian(TS ts,Mat Amat,PetscErrorCode (*func)(TS,PetscReal,Vec,Mat,void*),void *ctx)
1216814e85d6SHong Zhang {
1217814e85d6SHong Zhang   PetscErrorCode ierr;
1218814e85d6SHong Zhang 
1219814e85d6SHong Zhang   PetscFunctionBegin;
1220814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1221814e85d6SHong Zhang   PetscValidHeaderSpecific(Amat,MAT_CLASSID,2);
1222814e85d6SHong Zhang 
1223814e85d6SHong Zhang   ts->rhsjacobianp    = func;
1224814e85d6SHong Zhang   ts->rhsjacobianpctx = ctx;
1225814e85d6SHong Zhang   if(Amat) {
1226814e85d6SHong Zhang     ierr = PetscObjectReference((PetscObject)Amat);CHKERRQ(ierr);
1227814e85d6SHong Zhang     ierr = MatDestroy(&ts->Jacp);CHKERRQ(ierr);
1228814e85d6SHong Zhang     ts->Jacp = Amat;
1229814e85d6SHong Zhang   }
1230814e85d6SHong Zhang   PetscFunctionReturn(0);
1231814e85d6SHong Zhang }
1232814e85d6SHong Zhang 
1233814e85d6SHong Zhang /*@C
1234814e85d6SHong Zhang   TSAdjointComputeRHSJacobian - Deprecated, use TSComputeRHSJacobianP()
1235814e85d6SHong Zhang 
1236814e85d6SHong Zhang   Level: deprecated
1237814e85d6SHong Zhang 
1238814e85d6SHong Zhang @*/
1239c9ad9de0SHong Zhang PetscErrorCode TSAdjointComputeRHSJacobian(TS ts,PetscReal t,Vec U,Mat Amat)
1240814e85d6SHong Zhang {
1241814e85d6SHong Zhang   PetscErrorCode ierr;
1242814e85d6SHong Zhang 
1243814e85d6SHong Zhang   PetscFunctionBegin;
1244814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1245c9ad9de0SHong Zhang   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
1246814e85d6SHong Zhang   PetscValidPointer(Amat,4);
1247814e85d6SHong Zhang 
1248814e85d6SHong Zhang   PetscStackPush("TS user JacobianP function for sensitivity analysis");
1249c9ad9de0SHong Zhang   ierr = (*ts->rhsjacobianp)(ts,t,U,Amat,ts->rhsjacobianpctx); CHKERRQ(ierr);
1250814e85d6SHong Zhang   PetscStackPop;
1251814e85d6SHong Zhang   PetscFunctionReturn(0);
1252814e85d6SHong Zhang }
1253814e85d6SHong Zhang 
1254814e85d6SHong Zhang /*@
1255c9ad9de0SHong Zhang   TSAdjointComputeDRDYFunction - Deprecated, use TSComputeDRDUFunction()
1256814e85d6SHong Zhang 
1257814e85d6SHong Zhang   Level: deprecated
1258814e85d6SHong Zhang 
1259814e85d6SHong Zhang @*/
1260c9ad9de0SHong Zhang PetscErrorCode TSAdjointComputeDRDYFunction(TS ts,PetscReal t,Vec U,Vec *DRDU)
1261814e85d6SHong Zhang {
1262814e85d6SHong Zhang   PetscErrorCode ierr;
1263814e85d6SHong Zhang 
1264814e85d6SHong Zhang   PetscFunctionBegin;
1265814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1266c9ad9de0SHong Zhang   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
1267814e85d6SHong Zhang 
1268814e85d6SHong Zhang   PetscStackPush("TS user DRDY function for sensitivity analysis");
1269c9ad9de0SHong Zhang   ierr = (*ts->drdufunction)(ts,t,U,DRDU,ts->costintegrandctx); CHKERRQ(ierr);
1270814e85d6SHong Zhang   PetscStackPop;
1271814e85d6SHong Zhang   PetscFunctionReturn(0);
1272814e85d6SHong Zhang }
1273814e85d6SHong Zhang 
1274814e85d6SHong Zhang /*@
1275814e85d6SHong Zhang   TSAdjointComputeDRDPFunction - Deprecated, use TSComputeDRDPFunction()
1276814e85d6SHong Zhang 
1277814e85d6SHong Zhang   Level: deprecated
1278814e85d6SHong Zhang 
1279814e85d6SHong Zhang @*/
1280c9ad9de0SHong Zhang PetscErrorCode TSAdjointComputeDRDPFunction(TS ts,PetscReal t,Vec U,Vec *DRDP)
1281814e85d6SHong Zhang {
1282814e85d6SHong Zhang   PetscErrorCode ierr;
1283814e85d6SHong Zhang 
1284814e85d6SHong Zhang   PetscFunctionBegin;
1285814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1286c9ad9de0SHong Zhang   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
1287814e85d6SHong Zhang 
1288814e85d6SHong Zhang   PetscStackPush("TS user DRDP function for sensitivity analysis");
1289c9ad9de0SHong Zhang   ierr = (*ts->drdpfunction)(ts,t,U,DRDP,ts->costintegrandctx); CHKERRQ(ierr);
1290814e85d6SHong Zhang   PetscStackPop;
1291814e85d6SHong Zhang   PetscFunctionReturn(0);
1292814e85d6SHong Zhang }
1293814e85d6SHong Zhang 
1294814e85d6SHong Zhang /*@C
1295814e85d6SHong Zhang    TSAdjointMonitorSensi - monitors the first lambda sensitivity
1296814e85d6SHong Zhang 
1297814e85d6SHong Zhang    Level: intermediate
1298814e85d6SHong Zhang 
1299814e85d6SHong Zhang .keywords: TS, set, monitor
1300814e85d6SHong Zhang 
1301814e85d6SHong Zhang .seealso: TSAdjointMonitorSet()
1302814e85d6SHong Zhang @*/
1303814e85d6SHong Zhang PetscErrorCode TSAdjointMonitorSensi(TS ts,PetscInt step,PetscReal ptime,Vec v,PetscInt numcost,Vec *lambda,Vec *mu,PetscViewerAndFormat *vf)
1304814e85d6SHong Zhang {
1305814e85d6SHong Zhang   PetscErrorCode ierr;
1306814e85d6SHong Zhang   PetscViewer    viewer = vf->viewer;
1307814e85d6SHong Zhang 
1308814e85d6SHong Zhang   PetscFunctionBegin;
1309814e85d6SHong Zhang   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4);
1310814e85d6SHong Zhang   ierr = PetscViewerPushFormat(viewer,vf->format);CHKERRQ(ierr);
1311814e85d6SHong Zhang   ierr = VecView(lambda[0],viewer);CHKERRQ(ierr);
1312814e85d6SHong Zhang   ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
1313814e85d6SHong Zhang   PetscFunctionReturn(0);
1314814e85d6SHong Zhang }
1315814e85d6SHong Zhang 
1316814e85d6SHong Zhang /*@C
1317814e85d6SHong Zhang    TSAdjointMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user
1318814e85d6SHong Zhang 
1319814e85d6SHong Zhang    Collective on TS
1320814e85d6SHong Zhang 
1321814e85d6SHong Zhang    Input Parameters:
1322814e85d6SHong Zhang +  ts - TS object you wish to monitor
1323814e85d6SHong Zhang .  name - the monitor type one is seeking
1324814e85d6SHong Zhang .  help - message indicating what monitoring is done
1325814e85d6SHong Zhang .  manual - manual page for the monitor
1326814e85d6SHong Zhang .  monitor - the monitor function
1327814e85d6SHong Zhang -  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
1328814e85d6SHong Zhang 
1329814e85d6SHong Zhang    Level: developer
1330814e85d6SHong Zhang 
1331814e85d6SHong Zhang .seealso: PetscOptionsGetViewer(), PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(),
1332814e85d6SHong Zhang           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
1333814e85d6SHong Zhang           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
1334814e85d6SHong Zhang           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
1335814e85d6SHong Zhang           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1336814e85d6SHong Zhang           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1337814e85d6SHong Zhang           PetscOptionsFList(), PetscOptionsEList()
1338814e85d6SHong Zhang @*/
1339814e85d6SHong Zhang 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*))
1340814e85d6SHong Zhang {
1341814e85d6SHong Zhang   PetscErrorCode    ierr;
1342814e85d6SHong Zhang   PetscViewer       viewer;
1343814e85d6SHong Zhang   PetscViewerFormat format;
1344814e85d6SHong Zhang   PetscBool         flg;
1345814e85d6SHong Zhang 
1346814e85d6SHong Zhang   PetscFunctionBegin;
134716413a6aSBarry Smith   ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)ts),((PetscObject) ts)->options,((PetscObject)ts)->prefix,name,&viewer,&format,&flg);CHKERRQ(ierr);
1348814e85d6SHong Zhang   if (flg) {
1349814e85d6SHong Zhang     PetscViewerAndFormat *vf;
1350814e85d6SHong Zhang     ierr = PetscViewerAndFormatCreate(viewer,format,&vf);CHKERRQ(ierr);
1351814e85d6SHong Zhang     ierr = PetscObjectDereference((PetscObject)viewer);CHKERRQ(ierr);
1352814e85d6SHong Zhang     if (monitorsetup) {
1353814e85d6SHong Zhang       ierr = (*monitorsetup)(ts,vf);CHKERRQ(ierr);
1354814e85d6SHong Zhang     }
1355814e85d6SHong Zhang     ierr = TSAdjointMonitorSet(ts,(PetscErrorCode (*)(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*))monitor,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy);CHKERRQ(ierr);
1356814e85d6SHong Zhang   }
1357814e85d6SHong Zhang   PetscFunctionReturn(0);
1358814e85d6SHong Zhang }
1359814e85d6SHong Zhang 
1360814e85d6SHong Zhang /*@C
1361814e85d6SHong Zhang    TSAdjointMonitorSet - Sets an ADDITIONAL function that is to be used at every
1362814e85d6SHong Zhang    timestep to display the iteration's  progress.
1363814e85d6SHong Zhang 
1364814e85d6SHong Zhang    Logically Collective on TS
1365814e85d6SHong Zhang 
1366814e85d6SHong Zhang    Input Parameters:
1367814e85d6SHong Zhang +  ts - the TS context obtained from TSCreate()
1368814e85d6SHong Zhang .  adjointmonitor - monitoring routine
1369814e85d6SHong Zhang .  adjointmctx - [optional] user-defined context for private data for the
1370814e85d6SHong Zhang              monitor routine (use NULL if no context is desired)
1371814e85d6SHong Zhang -  adjointmonitordestroy - [optional] routine that frees monitor context
1372814e85d6SHong Zhang           (may be NULL)
1373814e85d6SHong Zhang 
1374814e85d6SHong Zhang    Calling sequence of monitor:
1375814e85d6SHong Zhang $    int adjointmonitor(TS ts,PetscInt steps,PetscReal time,Vec u,PetscInt numcost,Vec *lambda, Vec *mu,void *adjointmctx)
1376814e85d6SHong Zhang 
1377814e85d6SHong Zhang +    ts - the TS context
1378814e85d6SHong Zhang .    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
1379814e85d6SHong Zhang                                been interpolated to)
1380814e85d6SHong Zhang .    time - current time
1381814e85d6SHong Zhang .    u - current iterate
1382814e85d6SHong Zhang .    numcost - number of cost functionos
1383814e85d6SHong Zhang .    lambda - sensitivities to initial conditions
1384814e85d6SHong Zhang .    mu - sensitivities to parameters
1385814e85d6SHong Zhang -    adjointmctx - [optional] adjoint monitoring context
1386814e85d6SHong Zhang 
1387814e85d6SHong Zhang    Notes:
1388814e85d6SHong Zhang    This routine adds an additional monitor to the list of monitors that
1389814e85d6SHong Zhang    already has been loaded.
1390814e85d6SHong Zhang 
1391814e85d6SHong Zhang    Fortran Notes:
1392814e85d6SHong Zhang     Only a single monitor function can be set for each TS object
1393814e85d6SHong Zhang 
1394814e85d6SHong Zhang    Level: intermediate
1395814e85d6SHong Zhang 
1396814e85d6SHong Zhang .keywords: TS, timestep, set, adjoint, monitor
1397814e85d6SHong Zhang 
1398814e85d6SHong Zhang .seealso: TSAdjointMonitorCancel()
1399814e85d6SHong Zhang @*/
1400814e85d6SHong Zhang PetscErrorCode TSAdjointMonitorSet(TS ts,PetscErrorCode (*adjointmonitor)(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*),void *adjointmctx,PetscErrorCode (*adjointmdestroy)(void**))
1401814e85d6SHong Zhang {
1402814e85d6SHong Zhang   PetscErrorCode ierr;
1403814e85d6SHong Zhang   PetscInt       i;
1404814e85d6SHong Zhang   PetscBool      identical;
1405814e85d6SHong Zhang 
1406814e85d6SHong Zhang   PetscFunctionBegin;
1407814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1408814e85d6SHong Zhang   for (i=0; i<ts->numbermonitors;i++) {
1409814e85d6SHong Zhang     ierr = PetscMonitorCompare((PetscErrorCode (*)(void))adjointmonitor,adjointmctx,adjointmdestroy,(PetscErrorCode (*)(void))ts->adjointmonitor[i],ts->adjointmonitorcontext[i],ts->adjointmonitordestroy[i],&identical);CHKERRQ(ierr);
1410814e85d6SHong Zhang     if (identical) PetscFunctionReturn(0);
1411814e85d6SHong Zhang   }
1412814e85d6SHong Zhang   if (ts->numberadjointmonitors >= MAXTSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many adjoint monitors set");
1413814e85d6SHong Zhang   ts->adjointmonitor[ts->numberadjointmonitors]          = adjointmonitor;
1414814e85d6SHong Zhang   ts->adjointmonitordestroy[ts->numberadjointmonitors]   = adjointmdestroy;
1415814e85d6SHong Zhang   ts->adjointmonitorcontext[ts->numberadjointmonitors++] = (void*)adjointmctx;
1416814e85d6SHong Zhang   PetscFunctionReturn(0);
1417814e85d6SHong Zhang }
1418814e85d6SHong Zhang 
1419814e85d6SHong Zhang /*@C
1420814e85d6SHong Zhang    TSAdjointMonitorCancel - Clears all the adjoint monitors that have been set on a time-step object.
1421814e85d6SHong Zhang 
1422814e85d6SHong Zhang    Logically Collective on TS
1423814e85d6SHong Zhang 
1424814e85d6SHong Zhang    Input Parameters:
1425814e85d6SHong Zhang .  ts - the TS context obtained from TSCreate()
1426814e85d6SHong Zhang 
1427814e85d6SHong Zhang    Notes:
1428814e85d6SHong Zhang    There is no way to remove a single, specific monitor.
1429814e85d6SHong Zhang 
1430814e85d6SHong Zhang    Level: intermediate
1431814e85d6SHong Zhang 
1432814e85d6SHong Zhang .keywords: TS, timestep, set, adjoint, monitor
1433814e85d6SHong Zhang 
1434814e85d6SHong Zhang .seealso: TSAdjointMonitorSet()
1435814e85d6SHong Zhang @*/
1436814e85d6SHong Zhang PetscErrorCode TSAdjointMonitorCancel(TS ts)
1437814e85d6SHong Zhang {
1438814e85d6SHong Zhang   PetscErrorCode ierr;
1439814e85d6SHong Zhang   PetscInt       i;
1440814e85d6SHong Zhang 
1441814e85d6SHong Zhang   PetscFunctionBegin;
1442814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1443814e85d6SHong Zhang   for (i=0; i<ts->numberadjointmonitors; i++) {
1444814e85d6SHong Zhang     if (ts->adjointmonitordestroy[i]) {
1445814e85d6SHong Zhang       ierr = (*ts->adjointmonitordestroy[i])(&ts->adjointmonitorcontext[i]);CHKERRQ(ierr);
1446814e85d6SHong Zhang     }
1447814e85d6SHong Zhang   }
1448814e85d6SHong Zhang   ts->numberadjointmonitors = 0;
1449814e85d6SHong Zhang   PetscFunctionReturn(0);
1450814e85d6SHong Zhang }
1451814e85d6SHong Zhang 
1452814e85d6SHong Zhang /*@C
1453814e85d6SHong Zhang    TSAdjointMonitorDefault - the default monitor of adjoint computations
1454814e85d6SHong Zhang 
1455814e85d6SHong Zhang    Level: intermediate
1456814e85d6SHong Zhang 
1457814e85d6SHong Zhang .keywords: TS, set, monitor
1458814e85d6SHong Zhang 
1459814e85d6SHong Zhang .seealso: TSAdjointMonitorSet()
1460814e85d6SHong Zhang @*/
1461814e85d6SHong Zhang PetscErrorCode TSAdjointMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,PetscInt numcost,Vec *lambda,Vec *mu,PetscViewerAndFormat *vf)
1462814e85d6SHong Zhang {
1463814e85d6SHong Zhang   PetscErrorCode ierr;
1464814e85d6SHong Zhang   PetscViewer    viewer = vf->viewer;
1465814e85d6SHong Zhang 
1466814e85d6SHong Zhang   PetscFunctionBegin;
1467814e85d6SHong Zhang   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4);
1468814e85d6SHong Zhang   ierr = PetscViewerPushFormat(viewer,vf->format);CHKERRQ(ierr);
1469814e85d6SHong Zhang   ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
1470814e85d6SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"%D TS dt %g time %g%s",step,(double)ts->time_step,(double)ptime,ts->steprollback ? " (r)\n" : "\n");CHKERRQ(ierr);
1471814e85d6SHong Zhang   ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
1472814e85d6SHong Zhang   ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
1473814e85d6SHong Zhang   PetscFunctionReturn(0);
1474814e85d6SHong Zhang }
1475814e85d6SHong Zhang 
1476814e85d6SHong Zhang /*@C
1477814e85d6SHong Zhang    TSAdjointMonitorDrawSensi - Monitors progress of the adjoint TS solvers by calling
1478814e85d6SHong Zhang    VecView() for the sensitivities to initial states at each timestep
1479814e85d6SHong Zhang 
1480814e85d6SHong Zhang    Collective on TS
1481814e85d6SHong Zhang 
1482814e85d6SHong Zhang    Input Parameters:
1483814e85d6SHong Zhang +  ts - the TS context
1484814e85d6SHong Zhang .  step - current time-step
1485814e85d6SHong Zhang .  ptime - current time
1486814e85d6SHong Zhang .  u - current state
1487814e85d6SHong Zhang .  numcost - number of cost functions
1488814e85d6SHong Zhang .  lambda - sensitivities to initial conditions
1489814e85d6SHong Zhang .  mu - sensitivities to parameters
1490814e85d6SHong Zhang -  dummy - either a viewer or NULL
1491814e85d6SHong Zhang 
1492814e85d6SHong Zhang    Level: intermediate
1493814e85d6SHong Zhang 
1494814e85d6SHong Zhang .keywords: TS,  vector, adjoint, monitor, view
1495814e85d6SHong Zhang 
1496814e85d6SHong Zhang .seealso: TSAdjointMonitorSet(), TSAdjointMonitorDefault(), VecView()
1497814e85d6SHong Zhang @*/
1498814e85d6SHong Zhang PetscErrorCode TSAdjointMonitorDrawSensi(TS ts,PetscInt step,PetscReal ptime,Vec u,PetscInt numcost,Vec *lambda,Vec *mu,void *dummy)
1499814e85d6SHong Zhang {
1500814e85d6SHong Zhang   PetscErrorCode   ierr;
1501814e85d6SHong Zhang   TSMonitorDrawCtx ictx = (TSMonitorDrawCtx)dummy;
1502814e85d6SHong Zhang   PetscDraw        draw;
1503814e85d6SHong Zhang   PetscReal        xl,yl,xr,yr,h;
1504814e85d6SHong Zhang   char             time[32];
1505814e85d6SHong Zhang 
1506814e85d6SHong Zhang   PetscFunctionBegin;
1507814e85d6SHong Zhang   if (!(((ictx->howoften > 0) && (!(step % ictx->howoften))) || ((ictx->howoften == -1) && ts->reason))) PetscFunctionReturn(0);
1508814e85d6SHong Zhang 
1509814e85d6SHong Zhang   ierr = VecView(lambda[0],ictx->viewer);CHKERRQ(ierr);
1510814e85d6SHong Zhang   ierr = PetscViewerDrawGetDraw(ictx->viewer,0,&draw);CHKERRQ(ierr);
1511814e85d6SHong Zhang   ierr = PetscSNPrintf(time,32,"Timestep %d Time %g",(int)step,(double)ptime);CHKERRQ(ierr);
1512814e85d6SHong Zhang   ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
1513814e85d6SHong Zhang   h    = yl + .95*(yr - yl);
1514814e85d6SHong Zhang   ierr = PetscDrawStringCentered(draw,.5*(xl+xr),h,PETSC_DRAW_BLACK,time);CHKERRQ(ierr);
1515814e85d6SHong Zhang   ierr = PetscDrawFlush(draw);CHKERRQ(ierr);
1516814e85d6SHong Zhang   PetscFunctionReturn(0);
1517814e85d6SHong Zhang }
1518814e85d6SHong Zhang 
1519814e85d6SHong Zhang /*
1520814e85d6SHong Zhang    TSAdjointSetFromOptions - Sets various TSAdjoint parameters from user options.
1521814e85d6SHong Zhang 
1522814e85d6SHong Zhang    Collective on TSAdjoint
1523814e85d6SHong Zhang 
1524814e85d6SHong Zhang    Input Parameter:
1525814e85d6SHong Zhang .  ts - the TS context
1526814e85d6SHong Zhang 
1527814e85d6SHong Zhang    Options Database Keys:
1528814e85d6SHong Zhang +  -ts_adjoint_solve <yes,no> After solving the ODE/DAE solve the adjoint problem (requires -ts_save_trajectory)
1529814e85d6SHong Zhang .  -ts_adjoint_monitor - print information at each adjoint time step
1530814e85d6SHong Zhang -  -ts_adjoint_monitor_draw_sensi - monitor the sensitivity of the first cost function wrt initial conditions (lambda[0]) graphically
1531814e85d6SHong Zhang 
1532814e85d6SHong Zhang    Level: developer
1533814e85d6SHong Zhang 
1534814e85d6SHong Zhang    Notes:
1535814e85d6SHong Zhang     This is not normally called directly by users
1536814e85d6SHong Zhang 
1537814e85d6SHong Zhang .keywords: TS, trajectory, timestep, set, options, database
1538814e85d6SHong Zhang 
1539814e85d6SHong Zhang .seealso: TSSetSaveTrajectory(), TSTrajectorySetUp()
1540814e85d6SHong Zhang */
1541814e85d6SHong Zhang PetscErrorCode TSAdjointSetFromOptions(PetscOptionItems *PetscOptionsObject,TS ts)
1542814e85d6SHong Zhang {
1543814e85d6SHong Zhang   PetscBool      tflg,opt;
1544814e85d6SHong Zhang   PetscErrorCode ierr;
1545814e85d6SHong Zhang 
1546814e85d6SHong Zhang   PetscFunctionBegin;
1547814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,2);
1548814e85d6SHong Zhang   ierr = PetscOptionsHead(PetscOptionsObject,"TS Adjoint options");CHKERRQ(ierr);
1549814e85d6SHong Zhang   tflg = ts->adjoint_solve ? PETSC_TRUE : PETSC_FALSE;
1550814e85d6SHong Zhang   ierr = PetscOptionsBool("-ts_adjoint_solve","Solve the adjoint problem immediately after solving the forward problem","",tflg,&tflg,&opt);CHKERRQ(ierr);
1551814e85d6SHong Zhang   if (opt) {
1552814e85d6SHong Zhang     ierr = TSSetSaveTrajectory(ts);CHKERRQ(ierr);
1553814e85d6SHong Zhang     ts->adjoint_solve = tflg;
1554814e85d6SHong Zhang   }
1555814e85d6SHong Zhang   ierr = TSAdjointMonitorSetFromOptions(ts,"-ts_adjoint_monitor","Monitor adjoint timestep size","TSAdjointMonitorDefault",TSAdjointMonitorDefault,NULL);CHKERRQ(ierr);
1556814e85d6SHong Zhang   ierr = TSAdjointMonitorSetFromOptions(ts,"-ts_adjoint_monitor_sensi","Monitor sensitivity in the adjoint computation","TSAdjointMonitorSensi",TSAdjointMonitorSensi,NULL);CHKERRQ(ierr);
1557814e85d6SHong Zhang   opt  = PETSC_FALSE;
1558814e85d6SHong Zhang   ierr = PetscOptionsName("-ts_adjoint_monitor_draw_sensi","Monitor adjoint sensitivities (lambda only) graphically","TSAdjointMonitorDrawSensi",&opt);CHKERRQ(ierr);
1559814e85d6SHong Zhang   if (opt) {
1560814e85d6SHong Zhang     TSMonitorDrawCtx ctx;
1561814e85d6SHong Zhang     PetscInt         howoften = 1;
1562814e85d6SHong Zhang 
1563814e85d6SHong Zhang     ierr = PetscOptionsInt("-ts_adjoint_monitor_draw_sensi","Monitor adjoint sensitivities (lambda only) graphically","TSAdjointMonitorDrawSensi",howoften,&howoften,NULL);CHKERRQ(ierr);
1564814e85d6SHong Zhang     ierr = TSMonitorDrawCtxCreate(PetscObjectComm((PetscObject)ts),0,0,PETSC_DECIDE,PETSC_DECIDE,300,300,howoften,&ctx);CHKERRQ(ierr);
1565814e85d6SHong Zhang     ierr = TSAdjointMonitorSet(ts,TSAdjointMonitorDrawSensi,ctx,(PetscErrorCode (*)(void**))TSMonitorDrawCtxDestroy);CHKERRQ(ierr);
1566814e85d6SHong Zhang   }
1567814e85d6SHong Zhang   PetscFunctionReturn(0);
1568814e85d6SHong Zhang }
1569814e85d6SHong Zhang 
1570814e85d6SHong Zhang /*@
1571814e85d6SHong Zhang    TSAdjointStep - Steps one time step backward in the adjoint run
1572814e85d6SHong Zhang 
1573814e85d6SHong Zhang    Collective on TS
1574814e85d6SHong Zhang 
1575814e85d6SHong Zhang    Input Parameter:
1576814e85d6SHong Zhang .  ts - the TS context obtained from TSCreate()
1577814e85d6SHong Zhang 
1578814e85d6SHong Zhang    Level: intermediate
1579814e85d6SHong Zhang 
1580814e85d6SHong Zhang .keywords: TS, adjoint, step
1581814e85d6SHong Zhang 
1582814e85d6SHong Zhang .seealso: TSAdjointSetUp(), TSAdjointSolve()
1583814e85d6SHong Zhang @*/
1584814e85d6SHong Zhang PetscErrorCode TSAdjointStep(TS ts)
1585814e85d6SHong Zhang {
1586814e85d6SHong Zhang   DM               dm;
1587814e85d6SHong Zhang   PetscErrorCode   ierr;
1588814e85d6SHong Zhang 
1589814e85d6SHong Zhang   PetscFunctionBegin;
1590814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1591814e85d6SHong Zhang   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
1592814e85d6SHong Zhang   ierr = TSAdjointSetUp(ts);CHKERRQ(ierr);
1593814e85d6SHong Zhang 
1594814e85d6SHong Zhang   ts->reason = TS_CONVERGED_ITERATING;
1595814e85d6SHong Zhang   ts->ptime_prev = ts->ptime;
1596814e85d6SHong Zhang   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);
1597814e85d6SHong Zhang   ierr = PetscLogEventBegin(TS_AdjointStep,ts,0,0,0);CHKERRQ(ierr);
1598814e85d6SHong Zhang   ierr = (*ts->ops->adjointstep)(ts);CHKERRQ(ierr);
1599814e85d6SHong Zhang   ierr = PetscLogEventEnd(TS_AdjointStep,ts,0,0,0);CHKERRQ(ierr);
1600814e85d6SHong Zhang   ts->adjoint_steps++; ts->steps--;
1601814e85d6SHong Zhang 
1602814e85d6SHong Zhang   if (ts->reason < 0) {
1603814e85d6SHong Zhang     if (ts->errorifstepfailed) {
160405c9950eSHong Zhang       if (ts->reason == TSADJOINT_DIVERGED_LINEAR_SOLVE) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_NOT_CONVERGED,"TSAdjointStep has failed due to %s",TSConvergedReasons[ts->reason]);
160505c9950eSHong Zhang       else SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_NOT_CONVERGED,"TSAdjointStep has failed due to %s",TSConvergedReasons[ts->reason]);
1606814e85d6SHong Zhang     }
1607814e85d6SHong Zhang   } else if (!ts->reason) {
1608814e85d6SHong Zhang     if (ts->adjoint_steps >= ts->adjoint_max_steps) ts->reason = TS_CONVERGED_ITS;
1609814e85d6SHong Zhang   }
1610814e85d6SHong Zhang   PetscFunctionReturn(0);
1611814e85d6SHong Zhang }
1612814e85d6SHong Zhang 
1613814e85d6SHong Zhang /*@
1614814e85d6SHong Zhang    TSAdjointSolve - Solves the discrete ajoint problem for an ODE/DAE
1615814e85d6SHong Zhang 
1616814e85d6SHong Zhang    Collective on TS
1617814e85d6SHong Zhang 
1618814e85d6SHong Zhang    Input Parameter:
1619814e85d6SHong Zhang .  ts - the TS context obtained from TSCreate()
1620814e85d6SHong Zhang 
1621814e85d6SHong Zhang    Options Database:
1622814e85d6SHong Zhang . -ts_adjoint_view_solution <viewerinfo> - views the first gradient with respect to the initial values
1623814e85d6SHong Zhang 
1624814e85d6SHong Zhang    Level: intermediate
1625814e85d6SHong Zhang 
1626814e85d6SHong Zhang    Notes:
1627814e85d6SHong Zhang    This must be called after a call to TSSolve() that solves the forward problem
1628814e85d6SHong Zhang 
1629814e85d6SHong Zhang    By default this will integrate back to the initial time, one can use TSAdjointSetSteps() to step back to a later time
1630814e85d6SHong Zhang 
1631814e85d6SHong Zhang .keywords: TS, timestep, solve
1632814e85d6SHong Zhang 
1633814e85d6SHong Zhang .seealso: TSCreate(), TSSetCostGradients(), TSSetSolution(), TSAdjointStep()
1634814e85d6SHong Zhang @*/
1635814e85d6SHong Zhang PetscErrorCode TSAdjointSolve(TS ts)
1636814e85d6SHong Zhang {
1637814e85d6SHong Zhang   PetscErrorCode    ierr;
1638814e85d6SHong Zhang 
1639814e85d6SHong Zhang   PetscFunctionBegin;
1640814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1641814e85d6SHong Zhang   ierr = TSAdjointSetUp(ts);CHKERRQ(ierr);
1642814e85d6SHong Zhang 
1643814e85d6SHong Zhang   /* reset time step and iteration counters */
1644814e85d6SHong Zhang   ts->adjoint_steps     = 0;
1645814e85d6SHong Zhang   ts->ksp_its           = 0;
1646814e85d6SHong Zhang   ts->snes_its          = 0;
1647814e85d6SHong Zhang   ts->num_snes_failures = 0;
1648814e85d6SHong Zhang   ts->reject            = 0;
1649814e85d6SHong Zhang   ts->reason            = TS_CONVERGED_ITERATING;
1650814e85d6SHong Zhang 
1651814e85d6SHong Zhang   if (!ts->adjoint_max_steps) ts->adjoint_max_steps = ts->steps;
1652814e85d6SHong Zhang   if (ts->adjoint_steps >= ts->adjoint_max_steps) ts->reason = TS_CONVERGED_ITS;
1653814e85d6SHong Zhang 
1654814e85d6SHong Zhang   while (!ts->reason) {
1655814e85d6SHong Zhang     ierr = TSTrajectoryGet(ts->trajectory,ts,ts->steps,&ts->ptime);CHKERRQ(ierr);
1656814e85d6SHong Zhang     ierr = TSAdjointMonitor(ts,ts->steps,ts->ptime,ts->vec_sol,ts->numcost,ts->vecs_sensi,ts->vecs_sensip);CHKERRQ(ierr);
1657814e85d6SHong Zhang     ierr = TSAdjointEventHandler(ts);CHKERRQ(ierr);
1658814e85d6SHong Zhang     ierr = TSAdjointStep(ts);CHKERRQ(ierr);
1659cd4cee2dSHong Zhang     if ((ts->vec_costintegral || ts->quadraturets) && !ts->costintegralfwd) {
1660814e85d6SHong Zhang       ierr = TSAdjointCostIntegral(ts);CHKERRQ(ierr);
1661814e85d6SHong Zhang     }
1662814e85d6SHong Zhang   }
1663814e85d6SHong Zhang   ierr = TSTrajectoryGet(ts->trajectory,ts,ts->steps,&ts->ptime);CHKERRQ(ierr);
1664814e85d6SHong Zhang   ierr = TSAdjointMonitor(ts,ts->steps,ts->ptime,ts->vec_sol,ts->numcost,ts->vecs_sensi,ts->vecs_sensip);CHKERRQ(ierr);
1665814e85d6SHong Zhang   ts->solvetime = ts->ptime;
1666814e85d6SHong Zhang   ierr = TSTrajectoryViewFromOptions(ts->trajectory,NULL,"-ts_trajectory_view");CHKERRQ(ierr);
1667814e85d6SHong Zhang   ierr = VecViewFromOptions(ts->vecs_sensi[0],(PetscObject) ts, "-ts_adjoint_view_solution");CHKERRQ(ierr);
1668814e85d6SHong Zhang   ts->adjoint_max_steps = 0;
1669814e85d6SHong Zhang   PetscFunctionReturn(0);
1670814e85d6SHong Zhang }
1671814e85d6SHong Zhang 
1672814e85d6SHong Zhang /*@C
1673814e85d6SHong Zhang    TSAdjointMonitor - Runs all user-provided adjoint monitor routines set using TSAdjointMonitorSet()
1674814e85d6SHong Zhang 
1675814e85d6SHong Zhang    Collective on TS
1676814e85d6SHong Zhang 
1677814e85d6SHong Zhang    Input Parameters:
1678814e85d6SHong Zhang +  ts - time stepping context obtained from TSCreate()
1679814e85d6SHong Zhang .  step - step number that has just completed
1680814e85d6SHong Zhang .  ptime - model time of the state
1681814e85d6SHong Zhang .  u - state at the current model time
1682814e85d6SHong Zhang .  numcost - number of cost functions (dimension of lambda  or mu)
1683814e85d6SHong Zhang .  lambda - vectors containing the gradients of the cost functions with respect to the ODE/DAE solution variables
1684814e85d6SHong Zhang -  mu - vectors containing the gradients of the cost functions with respect to the problem parameters
1685814e85d6SHong Zhang 
1686814e85d6SHong Zhang    Notes:
1687814e85d6SHong Zhang    TSAdjointMonitor() is typically used automatically within the time stepping implementations.
1688814e85d6SHong Zhang    Users would almost never call this routine directly.
1689814e85d6SHong Zhang 
1690814e85d6SHong Zhang    Level: developer
1691814e85d6SHong Zhang 
1692814e85d6SHong Zhang .keywords: TS, timestep
1693814e85d6SHong Zhang @*/
1694814e85d6SHong Zhang PetscErrorCode TSAdjointMonitor(TS ts,PetscInt step,PetscReal ptime,Vec u,PetscInt numcost,Vec *lambda, Vec *mu)
1695814e85d6SHong Zhang {
1696814e85d6SHong Zhang   PetscErrorCode ierr;
1697814e85d6SHong Zhang   PetscInt       i,n = ts->numberadjointmonitors;
1698814e85d6SHong Zhang 
1699814e85d6SHong Zhang   PetscFunctionBegin;
1700814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1701814e85d6SHong Zhang   PetscValidHeaderSpecific(u,VEC_CLASSID,4);
17028860a134SJunchao Zhang   ierr = VecLockReadPush(u);CHKERRQ(ierr);
1703814e85d6SHong Zhang   for (i=0; i<n; i++) {
1704814e85d6SHong Zhang     ierr = (*ts->adjointmonitor[i])(ts,step,ptime,u,numcost,lambda,mu,ts->adjointmonitorcontext[i]);CHKERRQ(ierr);
1705814e85d6SHong Zhang   }
17068860a134SJunchao Zhang   ierr = VecLockReadPop(u);CHKERRQ(ierr);
1707814e85d6SHong Zhang   PetscFunctionReturn(0);
1708814e85d6SHong Zhang }
1709814e85d6SHong Zhang 
1710814e85d6SHong Zhang /*@
1711814e85d6SHong Zhang  TSAdjointCostIntegral - Evaluate the cost integral in the adjoint run.
1712814e85d6SHong Zhang 
1713814e85d6SHong Zhang  Collective on TS
1714814e85d6SHong Zhang 
1715814e85d6SHong Zhang  Input Arguments:
1716814e85d6SHong Zhang  .  ts - time stepping context
1717814e85d6SHong Zhang 
1718814e85d6SHong Zhang  Level: advanced
1719814e85d6SHong Zhang 
1720814e85d6SHong Zhang  Notes:
1721814e85d6SHong Zhang  This function cannot be called until TSAdjointStep() has been completed.
1722814e85d6SHong Zhang 
1723814e85d6SHong Zhang  .seealso: TSAdjointSolve(), TSAdjointStep
1724814e85d6SHong Zhang  @*/
1725814e85d6SHong Zhang PetscErrorCode TSAdjointCostIntegral(TS ts)
1726814e85d6SHong Zhang {
1727814e85d6SHong Zhang     PetscErrorCode ierr;
1728814e85d6SHong Zhang     PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1729814e85d6SHong Zhang     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);
1730814e85d6SHong Zhang     ierr = (*ts->ops->adjointintegral)(ts);CHKERRQ(ierr);
1731814e85d6SHong Zhang     PetscFunctionReturn(0);
1732814e85d6SHong Zhang }
1733814e85d6SHong Zhang 
1734814e85d6SHong Zhang /* ------------------ Forward (tangent linear) sensitivity  ------------------*/
1735814e85d6SHong Zhang 
1736814e85d6SHong Zhang /*@
1737814e85d6SHong Zhang   TSForwardSetUp - Sets up the internal data structures for the later use
1738814e85d6SHong Zhang   of forward sensitivity analysis
1739814e85d6SHong Zhang 
1740814e85d6SHong Zhang   Collective on TS
1741814e85d6SHong Zhang 
1742814e85d6SHong Zhang   Input Parameter:
1743814e85d6SHong Zhang . ts - the TS context obtained from TSCreate()
1744814e85d6SHong Zhang 
1745814e85d6SHong Zhang   Level: advanced
1746814e85d6SHong Zhang 
1747814e85d6SHong Zhang .keywords: TS, forward sensitivity, setup
1748814e85d6SHong Zhang 
1749814e85d6SHong Zhang .seealso: TSCreate(), TSDestroy(), TSSetUp()
1750814e85d6SHong Zhang @*/
1751814e85d6SHong Zhang PetscErrorCode TSForwardSetUp(TS ts)
1752814e85d6SHong Zhang {
1753814e85d6SHong Zhang   PetscErrorCode ierr;
1754814e85d6SHong Zhang 
1755814e85d6SHong Zhang   PetscFunctionBegin;
1756814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1757814e85d6SHong Zhang   if (ts->forwardsetupcalled) PetscFunctionReturn(0);
175833f72c5dSHong Zhang   if (!ts->Jacp && ts->Jacprhs) ts->Jacp = ts->Jacprhs;
1759814e85d6SHong Zhang 
1760814e85d6SHong Zhang   if (ts->ops->forwardsetup) {
1761814e85d6SHong Zhang     ierr = (*ts->ops->forwardsetup)(ts);CHKERRQ(ierr);
1762814e85d6SHong Zhang   }
17637207e0fdSHong Zhang   ierr = VecDuplicate(ts->vec_sol,&ts->vec_sensip_col);CHKERRQ(ierr);
1764814e85d6SHong Zhang   ts->forwardsetupcalled = PETSC_TRUE;
1765814e85d6SHong Zhang   PetscFunctionReturn(0);
1766814e85d6SHong Zhang }
1767814e85d6SHong Zhang 
1768814e85d6SHong Zhang /*@
17697adebddeSHong Zhang   TSForwardReset - Reset the internal data structures used by forward sensitivity analysis
17707adebddeSHong Zhang 
17717adebddeSHong Zhang   Collective on TS
17727adebddeSHong Zhang 
17737adebddeSHong Zhang   Input Parameter:
17747adebddeSHong Zhang . ts - the TS context obtained from TSCreate()
17757adebddeSHong Zhang 
17767adebddeSHong Zhang   Level: advanced
17777adebddeSHong Zhang 
17787adebddeSHong Zhang .keywords: TS, forward sensitivity, reset
17797adebddeSHong Zhang 
17807adebddeSHong Zhang .seealso: TSCreate(), TSDestroy(), TSForwardSetUp()
17817adebddeSHong Zhang @*/
17827adebddeSHong Zhang PetscErrorCode TSForwardReset(TS ts)
17837adebddeSHong Zhang {
17847207e0fdSHong Zhang   TS             quadts = ts->quadraturets;
17857adebddeSHong Zhang   PetscErrorCode ierr;
17867adebddeSHong Zhang 
17877adebddeSHong Zhang   PetscFunctionBegin;
17887adebddeSHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
17897adebddeSHong Zhang   if (ts->ops->forwardreset) {
17907adebddeSHong Zhang     ierr = (*ts->ops->forwardreset)(ts);CHKERRQ(ierr);
17917adebddeSHong Zhang   }
17927adebddeSHong Zhang   ierr = MatDestroy(&ts->mat_sensip);CHKERRQ(ierr);
17937207e0fdSHong Zhang   if (quadts) {
17947207e0fdSHong Zhang     ierr = MatDestroy(&quadts->mat_sensip);CHKERRQ(ierr);
17957207e0fdSHong Zhang   }
17967207e0fdSHong Zhang   ierr = VecDestroy(&ts->vec_sensip_col);CHKERRQ(ierr);
17977adebddeSHong Zhang   ts->forward_solve      = PETSC_FALSE;
17987adebddeSHong Zhang   ts->forwardsetupcalled = PETSC_FALSE;
17997adebddeSHong Zhang   PetscFunctionReturn(0);
18007adebddeSHong Zhang }
18017adebddeSHong Zhang 
18027adebddeSHong Zhang /*@
1803814e85d6SHong Zhang   TSForwardSetIntegralGradients - Set the vectors holding forward sensitivities of the integral term.
1804814e85d6SHong Zhang 
1805814e85d6SHong Zhang   Input Parameter:
1806814e85d6SHong Zhang . ts- the TS context obtained from TSCreate()
1807814e85d6SHong Zhang . numfwdint- number of integrals
1808814e85d6SHong Zhang . vp = the vectors containing the gradients for each integral w.r.t. parameters
1809814e85d6SHong Zhang 
18107207e0fdSHong Zhang   Level: deprecated
1811814e85d6SHong Zhang 
1812814e85d6SHong Zhang .keywords: TS, forward sensitivity
1813814e85d6SHong Zhang 
1814814e85d6SHong Zhang .seealso: TSForwardGetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep()
1815814e85d6SHong Zhang @*/
1816814e85d6SHong Zhang PetscErrorCode TSForwardSetIntegralGradients(TS ts,PetscInt numfwdint,Vec *vp)
1817814e85d6SHong Zhang {
1818814e85d6SHong Zhang   PetscFunctionBegin;
1819814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1820814e85d6SHong Zhang   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()");
1821814e85d6SHong Zhang   if (!ts->numcost) ts->numcost = numfwdint;
1822814e85d6SHong Zhang 
1823814e85d6SHong Zhang   ts->vecs_integral_sensip = vp;
1824814e85d6SHong Zhang   PetscFunctionReturn(0);
1825814e85d6SHong Zhang }
1826814e85d6SHong Zhang 
1827814e85d6SHong Zhang /*@
1828814e85d6SHong Zhang   TSForwardGetIntegralGradients - Returns the forward sensitivities ofthe integral term.
1829814e85d6SHong Zhang 
1830814e85d6SHong Zhang   Input Parameter:
1831814e85d6SHong Zhang . ts- the TS context obtained from TSCreate()
1832814e85d6SHong Zhang 
1833814e85d6SHong Zhang   Output Parameter:
1834814e85d6SHong Zhang . vp = the vectors containing the gradients for each integral w.r.t. parameters
1835814e85d6SHong Zhang 
18367207e0fdSHong Zhang   Level: deprecated
1837814e85d6SHong Zhang 
1838814e85d6SHong Zhang .keywords: TS, forward sensitivity
1839814e85d6SHong Zhang 
1840814e85d6SHong Zhang .seealso: TSForwardSetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep()
1841814e85d6SHong Zhang @*/
1842814e85d6SHong Zhang PetscErrorCode TSForwardGetIntegralGradients(TS ts,PetscInt *numfwdint,Vec **vp)
1843814e85d6SHong Zhang {
1844814e85d6SHong Zhang   PetscFunctionBegin;
1845814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1846814e85d6SHong Zhang   PetscValidPointer(vp,3);
1847814e85d6SHong Zhang   if (numfwdint) *numfwdint = ts->numcost;
1848814e85d6SHong Zhang   if (vp) *vp = ts->vecs_integral_sensip;
1849814e85d6SHong Zhang   PetscFunctionReturn(0);
1850814e85d6SHong Zhang }
1851814e85d6SHong Zhang 
1852814e85d6SHong Zhang /*@
1853814e85d6SHong Zhang   TSForwardStep - Compute the forward sensitivity for one time step.
1854814e85d6SHong Zhang 
1855814e85d6SHong Zhang   Collective on TS
1856814e85d6SHong Zhang 
1857814e85d6SHong Zhang   Input Arguments:
1858814e85d6SHong Zhang . ts - time stepping context
1859814e85d6SHong Zhang 
1860814e85d6SHong Zhang   Level: advanced
1861814e85d6SHong Zhang 
1862814e85d6SHong Zhang   Notes:
1863814e85d6SHong Zhang   This function cannot be called until TSStep() has been completed.
1864814e85d6SHong Zhang 
1865814e85d6SHong Zhang .keywords: TS, forward sensitivity
1866814e85d6SHong Zhang 
1867814e85d6SHong Zhang .seealso: TSForwardSetSensitivities(), TSForwardGetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardSetUp()
1868814e85d6SHong Zhang @*/
1869814e85d6SHong Zhang PetscErrorCode TSForwardStep(TS ts)
1870814e85d6SHong Zhang {
1871814e85d6SHong Zhang   PetscErrorCode ierr;
1872814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1873814e85d6SHong Zhang   if (!ts->ops->forwardstep) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"%s does not provide forward sensitivity analysis",((PetscObject)ts)->type_name);
1874814e85d6SHong Zhang   ierr = PetscLogEventBegin(TS_ForwardStep,ts,0,0,0);CHKERRQ(ierr);
1875814e85d6SHong Zhang   ierr = (*ts->ops->forwardstep)(ts);CHKERRQ(ierr);
1876814e85d6SHong Zhang   ierr = PetscLogEventEnd(TS_ForwardStep,ts,0,0,0);CHKERRQ(ierr);
1877814e85d6SHong Zhang   PetscFunctionReturn(0);
1878814e85d6SHong Zhang }
1879814e85d6SHong Zhang 
1880814e85d6SHong Zhang /*@
1881814e85d6SHong Zhang   TSForwardSetSensitivities - Sets the initial value of the trajectory sensitivities of solution  w.r.t. the problem parameters and initial values.
1882814e85d6SHong Zhang 
1883814e85d6SHong Zhang   Logically Collective on TS and Vec
1884814e85d6SHong Zhang 
1885814e85d6SHong Zhang   Input Parameters:
1886814e85d6SHong Zhang + ts - the TS context obtained from TSCreate()
1887814e85d6SHong Zhang . nump - number of parameters
1888814e85d6SHong Zhang - Smat - sensitivities with respect to the parameters, the number of entries in these vectors is the same as the number of parameters
1889814e85d6SHong Zhang 
1890814e85d6SHong Zhang   Level: beginner
1891814e85d6SHong Zhang 
1892814e85d6SHong Zhang   Notes:
1893814e85d6SHong Zhang   Forward sensitivity is also called 'trajectory sensitivity' in some fields such as power systems.
1894814e85d6SHong Zhang   This function turns on a flag to trigger TSSolve() to compute forward sensitivities automatically.
1895814e85d6SHong Zhang   You must call this function before TSSolve().
1896814e85d6SHong Zhang   The entries in the sensitivity matrix must be correctly initialized with the values S = dy/dp|startingtime.
1897814e85d6SHong Zhang 
1898814e85d6SHong Zhang .keywords: TS, timestep, set, forward sensitivity, initial values
1899814e85d6SHong Zhang 
1900814e85d6SHong Zhang .seealso: TSForwardGetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep()
1901814e85d6SHong Zhang @*/
1902814e85d6SHong Zhang PetscErrorCode TSForwardSetSensitivities(TS ts,PetscInt nump,Mat Smat)
1903814e85d6SHong Zhang {
1904814e85d6SHong Zhang   PetscErrorCode ierr;
1905814e85d6SHong Zhang 
1906814e85d6SHong Zhang   PetscFunctionBegin;
1907814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1908814e85d6SHong Zhang   PetscValidHeaderSpecific(Smat,MAT_CLASSID,3);
1909814e85d6SHong Zhang   ts->forward_solve  = PETSC_TRUE;
1910b5b4017aSHong Zhang   if (nump == PETSC_DEFAULT) {
1911b5b4017aSHong Zhang     ierr = MatGetSize(Smat,NULL,&ts->num_parameters);CHKERRQ(ierr);
1912b5b4017aSHong Zhang   } else ts->num_parameters = nump;
1913814e85d6SHong Zhang   ierr = PetscObjectReference((PetscObject)Smat);CHKERRQ(ierr);
1914814e85d6SHong Zhang   ierr = MatDestroy(&ts->mat_sensip);CHKERRQ(ierr);
1915814e85d6SHong Zhang   ts->mat_sensip = Smat;
1916814e85d6SHong Zhang   PetscFunctionReturn(0);
1917814e85d6SHong Zhang }
1918814e85d6SHong Zhang 
1919814e85d6SHong Zhang /*@
1920814e85d6SHong Zhang   TSForwardGetSensitivities - Returns the trajectory sensitivities
1921814e85d6SHong Zhang 
1922814e85d6SHong Zhang   Not Collective, but Vec returned is parallel if TS is parallel
1923814e85d6SHong Zhang 
1924814e85d6SHong Zhang   Output Parameter:
1925814e85d6SHong Zhang + ts - the TS context obtained from TSCreate()
1926814e85d6SHong Zhang . nump - number of parameters
1927814e85d6SHong Zhang - Smat - sensitivities with respect to the parameters, the number of entries in these vectors is the same as the number of parameters
1928814e85d6SHong Zhang 
1929814e85d6SHong Zhang   Level: intermediate
1930814e85d6SHong Zhang 
1931814e85d6SHong Zhang .keywords: TS, forward sensitivity
1932814e85d6SHong Zhang 
1933814e85d6SHong Zhang .seealso: TSForwardSetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep()
1934814e85d6SHong Zhang @*/
1935814e85d6SHong Zhang PetscErrorCode TSForwardGetSensitivities(TS ts,PetscInt *nump,Mat *Smat)
1936814e85d6SHong Zhang {
1937814e85d6SHong Zhang   PetscFunctionBegin;
1938814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1939814e85d6SHong Zhang   if (nump) *nump = ts->num_parameters;
1940814e85d6SHong Zhang   if (Smat) *Smat = ts->mat_sensip;
1941814e85d6SHong Zhang   PetscFunctionReturn(0);
1942814e85d6SHong Zhang }
1943814e85d6SHong Zhang 
1944814e85d6SHong Zhang /*@
1945814e85d6SHong Zhang    TSForwardCostIntegral - Evaluate the cost integral in the forward run.
1946814e85d6SHong Zhang 
1947814e85d6SHong Zhang    Collective on TS
1948814e85d6SHong Zhang 
1949814e85d6SHong Zhang    Input Arguments:
1950814e85d6SHong Zhang .  ts - time stepping context
1951814e85d6SHong Zhang 
1952814e85d6SHong Zhang    Level: advanced
1953814e85d6SHong Zhang 
1954814e85d6SHong Zhang    Notes:
1955814e85d6SHong Zhang    This function cannot be called until TSStep() has been completed.
1956814e85d6SHong Zhang 
1957814e85d6SHong Zhang .seealso: TSSolve(), TSAdjointCostIntegral()
1958814e85d6SHong Zhang @*/
1959814e85d6SHong Zhang PetscErrorCode TSForwardCostIntegral(TS ts)
1960814e85d6SHong Zhang {
1961814e85d6SHong Zhang   PetscErrorCode ierr;
1962814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1963814e85d6SHong Zhang   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);
1964814e85d6SHong Zhang   ierr = (*ts->ops->forwardintegral)(ts);CHKERRQ(ierr);
1965814e85d6SHong Zhang   PetscFunctionReturn(0);
1966814e85d6SHong Zhang }
1967b5b4017aSHong Zhang 
1968b5b4017aSHong Zhang /*@
1969b5b4017aSHong Zhang   TSForwardSetInitialSensitivities - Set initial values for tangent linear sensitivities
1970b5b4017aSHong Zhang 
1971b5b4017aSHong Zhang   Collective on TS and Mat
1972b5b4017aSHong Zhang 
1973b5b4017aSHong Zhang   Input Parameter
1974b5b4017aSHong Zhang + ts - the TS context obtained from TSCreate()
1975b5b4017aSHong Zhang - didp - parametric sensitivities of the initial condition
1976b5b4017aSHong Zhang 
1977b5b4017aSHong Zhang   Level: intermediate
1978b5b4017aSHong Zhang 
1979b5b4017aSHong Zhang   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.
1980b5b4017aSHong Zhang 
1981b5b4017aSHong Zhang .seealso: TSForwardSetSensitivities()
1982b5b4017aSHong Zhang @*/
1983b5b4017aSHong Zhang PetscErrorCode TSForwardSetInitialSensitivities(TS ts,Mat didp)
1984b5b4017aSHong Zhang {
1985b5b4017aSHong Zhang   PetscErrorCode ierr;
1986b5b4017aSHong Zhang 
1987b5b4017aSHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1988b5b4017aSHong Zhang   PetscValidHeaderSpecific(didp,MAT_CLASSID,2);
1989b5b4017aSHong Zhang   if (!ts->mat_sensip) {
1990b5b4017aSHong Zhang     ierr = TSForwardSetSensitivities(ts,PETSC_DEFAULT,didp);CHKERRQ(ierr);
1991b5b4017aSHong Zhang   }
1992b5b4017aSHong Zhang   PetscFunctionReturn(0);
1993b5b4017aSHong Zhang }
19946affb6f8SHong Zhang 
19956affb6f8SHong Zhang /*@
19966affb6f8SHong Zhang    TSForwardGetStages - Get the number of stages and the tangent linear sensitivities at the intermediate stages
19976affb6f8SHong Zhang 
19986affb6f8SHong Zhang    Input Parameter:
19996affb6f8SHong Zhang .  ts - the TS context obtained from TSCreate()
20006affb6f8SHong Zhang 
20016affb6f8SHong Zhang    Output Parameters:
2002cd4cee2dSHong Zhang +  ns - number of stages
20036affb6f8SHong Zhang -  S - tangent linear sensitivities at the intermediate stages
20046affb6f8SHong Zhang 
20056affb6f8SHong Zhang    Level: advanced
20066affb6f8SHong Zhang 
20076affb6f8SHong Zhang .keywords: TS, second-order adjoint, forward sensitivity
20086affb6f8SHong Zhang @*/
20096affb6f8SHong Zhang PetscErrorCode TSForwardGetStages(TS ts,PetscInt *ns,Mat **S)
20106affb6f8SHong Zhang {
20116affb6f8SHong Zhang   PetscErrorCode ierr;
20126affb6f8SHong Zhang 
20136affb6f8SHong Zhang   PetscFunctionBegin;
20146affb6f8SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
20156affb6f8SHong Zhang 
20166affb6f8SHong Zhang   if (!ts->ops->getstages) *S=NULL;
20176affb6f8SHong Zhang   else {
20186affb6f8SHong Zhang     ierr = (*ts->ops->forwardgetstages)(ts,ns,S);CHKERRQ(ierr);
20196affb6f8SHong Zhang   }
20206affb6f8SHong Zhang   PetscFunctionReturn(0);
20216affb6f8SHong Zhang }
2022cd4cee2dSHong Zhang 
2023cd4cee2dSHong Zhang /*@
2024cd4cee2dSHong Zhang    TSCreateQuadratureTS - Create a sub-TS that evaluates integrals over time
2025cd4cee2dSHong Zhang 
2026cd4cee2dSHong Zhang    Input Parameter:
2027cd4cee2dSHong Zhang +  ts - the TS context obtained from TSCreate()
2028cd4cee2dSHong Zhang -  fwd - flag indicating whether to evaluate cost integral in the forward run or the adjoint run
2029cd4cee2dSHong Zhang 
2030cd4cee2dSHong Zhang    Output Parameters:
2031cd4cee2dSHong Zhang .  quadts - the child TS context
2032cd4cee2dSHong Zhang 
2033cd4cee2dSHong Zhang    Level: intermediate
2034cd4cee2dSHong Zhang 
2035cd4cee2dSHong Zhang .keywords: TS, quadrature evaluation
2036cd4cee2dSHong Zhang 
2037cd4cee2dSHong Zhang .seealso: TSGetQuadratureTS()
2038cd4cee2dSHong Zhang @*/
2039cd4cee2dSHong Zhang PetscErrorCode TSCreateQuadratureTS(TS ts,PetscBool fwd,TS *quadts)
2040cd4cee2dSHong Zhang {
2041cd4cee2dSHong Zhang   char prefix[128];
2042cd4cee2dSHong Zhang   PetscErrorCode ierr;
2043cd4cee2dSHong Zhang 
2044cd4cee2dSHong Zhang   PetscFunctionBegin;
2045cd4cee2dSHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2046cd4cee2dSHong Zhang   PetscValidPointer(quadts,2);
2047cd4cee2dSHong Zhang   ierr = TSDestroy(&ts->quadraturets);CHKERRQ(ierr);
2048cd4cee2dSHong Zhang   ierr = TSCreate(PetscObjectComm((PetscObject)ts),&ts->quadraturets);CHKERRQ(ierr);
2049cd4cee2dSHong Zhang   ierr = PetscObjectIncrementTabLevel((PetscObject)ts->quadraturets,(PetscObject)ts,1);CHKERRQ(ierr);
2050cd4cee2dSHong Zhang   ierr = PetscLogObjectParent((PetscObject)ts,(PetscObject)ts->quadraturets);CHKERRQ(ierr);
2051cd4cee2dSHong Zhang   ierr = PetscSNPrintf(prefix,sizeof(prefix),"%squad_",((PetscObject)ts)->prefix ? ((PetscObject)ts)->prefix : "");CHKERRQ(ierr);
2052cd4cee2dSHong Zhang   ierr = TSSetOptionsPrefix(ts->quadraturets,prefix);CHKERRQ(ierr);
2053cd4cee2dSHong Zhang   *quadts = ts->quadraturets;
2054cd4cee2dSHong Zhang 
2055cd4cee2dSHong Zhang   if (ts->numcost) {
2056cd4cee2dSHong Zhang     ierr = VecCreateSeq(PETSC_COMM_SELF,ts->numcost,&(*quadts)->vec_sol);CHKERRQ(ierr);
2057cd4cee2dSHong Zhang   } else {
2058cd4cee2dSHong Zhang     ierr = VecCreateSeq(PETSC_COMM_SELF,1,&(*quadts)->vec_sol);CHKERRQ(ierr);
2059cd4cee2dSHong Zhang   }
2060cd4cee2dSHong Zhang   ts->costintegralfwd = fwd;
2061cd4cee2dSHong Zhang   PetscFunctionReturn(0);
2062cd4cee2dSHong Zhang }
2063cd4cee2dSHong Zhang 
2064cd4cee2dSHong Zhang /*@
2065cd4cee2dSHong Zhang    TSGetQuadratureTS - Return the sub-TS that evaluates integrals over time
2066cd4cee2dSHong Zhang 
2067cd4cee2dSHong Zhang    Input Parameter:
2068cd4cee2dSHong Zhang .  ts - the TS context obtained from TSCreate()
2069cd4cee2dSHong Zhang 
2070cd4cee2dSHong Zhang    Output Parameters:
2071cd4cee2dSHong Zhang +  fwd - flag indicating whether to evaluate cost integral in the forward run or the adjoint run
2072cd4cee2dSHong Zhang -  quadts - the child TS context
2073cd4cee2dSHong Zhang 
2074cd4cee2dSHong Zhang    Level: intermediate
2075cd4cee2dSHong Zhang 
2076cd4cee2dSHong Zhang .keywords: TS, quadrature evaluation
2077cd4cee2dSHong Zhang 
2078cd4cee2dSHong Zhang .seealso: TSCreateQuadratureTS()
2079cd4cee2dSHong Zhang @*/
2080cd4cee2dSHong Zhang PetscErrorCode TSGetQuadratureTS(TS ts,PetscBool *fwd,TS *quadts)
2081cd4cee2dSHong Zhang {
2082cd4cee2dSHong Zhang   PetscFunctionBegin;
2083cd4cee2dSHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
2084cd4cee2dSHong Zhang   if (fwd) *fwd = ts->costintegralfwd;
2085cd4cee2dSHong Zhang   if (quadts) *quadts = ts->quadraturets;
2086cd4cee2dSHong Zhang   PetscFunctionReturn(0);
2087cd4cee2dSHong Zhang }
2088