xref: /petsc/src/ts/interface/sensitivity/tssen.c (revision ba0675f606ccd55dd0554eb357639a24b77dc41b)
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 
67f59fb53SHong Zhang /* #define TSADJOINT_STAGE */
77f59fb53SHong Zhang 
8814e85d6SHong Zhang /* ------------------------ Sensitivity Context ---------------------------*/
9814e85d6SHong Zhang 
10814e85d6SHong Zhang /*@C
11b5b4017aSHong 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.
12814e85d6SHong Zhang 
13814e85d6SHong Zhang   Logically Collective on TS
14814e85d6SHong Zhang 
15814e85d6SHong Zhang   Input Parameters:
16b5b4017aSHong Zhang + ts - TS context obtained from TSCreate()
17b5b4017aSHong Zhang . Amat - JacobianP matrix
18b5b4017aSHong Zhang . func - function
19b5b4017aSHong Zhang - ctx - [optional] user-defined function context
20814e85d6SHong Zhang 
21814e85d6SHong Zhang   Calling sequence of func:
22814e85d6SHong Zhang $ func (TS ts,PetscReal t,Vec y,Mat A,void *ctx);
23814e85d6SHong Zhang +   t - current timestep
24c9ad9de0SHong Zhang .   U - input vector (current ODE solution)
25814e85d6SHong Zhang .   A - output matrix
26814e85d6SHong Zhang -   ctx - [optional] user-defined function context
27814e85d6SHong Zhang 
28814e85d6SHong Zhang   Level: intermediate
29814e85d6SHong Zhang 
30814e85d6SHong Zhang   Notes:
31814e85d6SHong 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
32814e85d6SHong Zhang 
33cd4cee2dSHong Zhang .seealso: TSGetRHSJacobianP()
34814e85d6SHong Zhang @*/
35814e85d6SHong Zhang PetscErrorCode TSSetRHSJacobianP(TS ts,Mat Amat,PetscErrorCode (*func)(TS,PetscReal,Vec,Mat,void*),void *ctx)
36814e85d6SHong Zhang {
37814e85d6SHong Zhang   PetscErrorCode ierr;
38814e85d6SHong Zhang 
39814e85d6SHong Zhang   PetscFunctionBegin;
40814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
41814e85d6SHong Zhang   PetscValidHeaderSpecific(Amat,MAT_CLASSID,2);
42814e85d6SHong Zhang 
43814e85d6SHong Zhang   ts->rhsjacobianp    = func;
44814e85d6SHong Zhang   ts->rhsjacobianpctx = ctx;
45814e85d6SHong Zhang   if(Amat) {
46814e85d6SHong Zhang     ierr = PetscObjectReference((PetscObject)Amat);CHKERRQ(ierr);
4733f72c5dSHong Zhang     ierr = MatDestroy(&ts->Jacprhs);CHKERRQ(ierr);
4833f72c5dSHong Zhang     ts->Jacprhs = Amat;
49814e85d6SHong Zhang   }
50814e85d6SHong Zhang   PetscFunctionReturn(0);
51814e85d6SHong Zhang }
52814e85d6SHong Zhang 
53814e85d6SHong Zhang /*@C
54cd4cee2dSHong 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.
55cd4cee2dSHong Zhang 
56cd4cee2dSHong Zhang   Logically Collective on TS
57cd4cee2dSHong Zhang 
58cd4cee2dSHong Zhang   Input Parameters:
59cd4cee2dSHong Zhang . ts - TS context obtained from TSCreate()
60cd4cee2dSHong Zhang 
61cd4cee2dSHong Zhang   Output Parameters:
62cd4cee2dSHong Zhang + Amat - JacobianP matrix
63cd4cee2dSHong Zhang . func - function
64cd4cee2dSHong Zhang - ctx - [optional] user-defined function context
65cd4cee2dSHong Zhang 
66cd4cee2dSHong Zhang   Calling sequence of func:
67cd4cee2dSHong Zhang $ func (TS ts,PetscReal t,Vec y,Mat A,void *ctx);
68cd4cee2dSHong Zhang +   t - current timestep
69cd4cee2dSHong Zhang .   U - input vector (current ODE solution)
70cd4cee2dSHong Zhang .   A - output matrix
71cd4cee2dSHong Zhang -   ctx - [optional] user-defined function context
72cd4cee2dSHong Zhang 
73cd4cee2dSHong Zhang   Level: intermediate
74cd4cee2dSHong Zhang 
75cd4cee2dSHong Zhang   Notes:
76cd4cee2dSHong 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
77cd4cee2dSHong Zhang 
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 .seealso: TSSetRHSJacobianP()
100814e85d6SHong Zhang @*/
101c9ad9de0SHong Zhang PetscErrorCode TSComputeRHSJacobianP(TS ts,PetscReal t,Vec U,Mat Amat)
102814e85d6SHong Zhang {
103814e85d6SHong Zhang   PetscErrorCode ierr;
104814e85d6SHong Zhang 
105814e85d6SHong Zhang   PetscFunctionBegin;
10633f72c5dSHong Zhang   if (!Amat) PetscFunctionReturn(0);
107814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
108c9ad9de0SHong Zhang   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
109814e85d6SHong Zhang 
110814e85d6SHong Zhang   PetscStackPush("TS user JacobianP function for sensitivity analysis");
111c9ad9de0SHong Zhang   ierr = (*ts->rhsjacobianp)(ts,t,U,Amat,ts->rhsjacobianpctx);CHKERRQ(ierr);
112814e85d6SHong Zhang   PetscStackPop;
113814e85d6SHong Zhang   PetscFunctionReturn(0);
114814e85d6SHong Zhang }
115814e85d6SHong Zhang 
116814e85d6SHong Zhang /*@C
11733f72c5dSHong 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.
11833f72c5dSHong Zhang 
11933f72c5dSHong Zhang   Logically Collective on TS
12033f72c5dSHong Zhang 
12133f72c5dSHong Zhang   Input Parameters:
12233f72c5dSHong Zhang + ts - TS context obtained from TSCreate()
12333f72c5dSHong Zhang . Amat - JacobianP matrix
12433f72c5dSHong Zhang . func - function
12533f72c5dSHong Zhang - ctx - [optional] user-defined function context
12633f72c5dSHong Zhang 
12733f72c5dSHong Zhang   Calling sequence of func:
12833f72c5dSHong Zhang $ func (TS ts,PetscReal t,Vec y,Mat A,void *ctx);
12933f72c5dSHong Zhang +   t - current timestep
13033f72c5dSHong Zhang .   U - input vector (current ODE solution)
13133f72c5dSHong Zhang .   Udot - time derivative of state vector
13233f72c5dSHong Zhang .   shift - shift to apply, see note below
13333f72c5dSHong Zhang .   A - output matrix
13433f72c5dSHong Zhang -   ctx - [optional] user-defined function context
13533f72c5dSHong Zhang 
13633f72c5dSHong Zhang   Level: intermediate
13733f72c5dSHong Zhang 
13833f72c5dSHong Zhang   Notes:
13933f72c5dSHong 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
14033f72c5dSHong Zhang 
14133f72c5dSHong Zhang .seealso:
14233f72c5dSHong Zhang @*/
14333f72c5dSHong Zhang PetscErrorCode TSSetIJacobianP(TS ts,Mat Amat,PetscErrorCode (*func)(TS,PetscReal,Vec,Vec,PetscReal,Mat,void*),void *ctx)
14433f72c5dSHong Zhang {
14533f72c5dSHong Zhang   PetscErrorCode ierr;
14633f72c5dSHong Zhang 
14733f72c5dSHong Zhang   PetscFunctionBegin;
14833f72c5dSHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
14933f72c5dSHong Zhang   PetscValidHeaderSpecific(Amat,MAT_CLASSID,2);
15033f72c5dSHong Zhang 
15133f72c5dSHong Zhang   ts->ijacobianp    = func;
15233f72c5dSHong Zhang   ts->ijacobianpctx = ctx;
15333f72c5dSHong Zhang   if(Amat) {
15433f72c5dSHong Zhang     ierr = PetscObjectReference((PetscObject)Amat);CHKERRQ(ierr);
15533f72c5dSHong Zhang     ierr = MatDestroy(&ts->Jacp);CHKERRQ(ierr);
15633f72c5dSHong Zhang     ts->Jacp = Amat;
15733f72c5dSHong Zhang   }
15833f72c5dSHong Zhang   PetscFunctionReturn(0);
15933f72c5dSHong Zhang }
16033f72c5dSHong Zhang 
16133f72c5dSHong Zhang /*@C
16233f72c5dSHong Zhang   TSComputeIJacobianP - Runs the user-defined IJacobianP function.
16333f72c5dSHong Zhang 
16433f72c5dSHong Zhang   Collective on TS
16533f72c5dSHong Zhang 
16633f72c5dSHong Zhang   Input Parameters:
16733f72c5dSHong Zhang + ts - the TS context
16833f72c5dSHong Zhang . t - current timestep
16933f72c5dSHong Zhang . U - state vector
17033f72c5dSHong Zhang . Udot - time derivative of state vector
17133f72c5dSHong Zhang . shift - shift to apply, see note below
17233f72c5dSHong Zhang - imex - flag indicates if the method is IMEX so that the RHSJacobian should be kept separate
17333f72c5dSHong Zhang 
17433f72c5dSHong Zhang   Output Parameters:
17533f72c5dSHong Zhang . A - Jacobian matrix
17633f72c5dSHong Zhang 
17733f72c5dSHong Zhang   Level: developer
17833f72c5dSHong Zhang 
17933f72c5dSHong Zhang .seealso: TSSetIJacobianP()
18033f72c5dSHong Zhang @*/
18133f72c5dSHong Zhang PetscErrorCode TSComputeIJacobianP(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal shift,Mat Amat,PetscBool imex)
18233f72c5dSHong Zhang {
18333f72c5dSHong Zhang   PetscErrorCode ierr;
18433f72c5dSHong Zhang 
18533f72c5dSHong Zhang   PetscFunctionBegin;
18633f72c5dSHong Zhang   if (!Amat) PetscFunctionReturn(0);
18733f72c5dSHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
18833f72c5dSHong Zhang   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
18933f72c5dSHong Zhang   PetscValidHeaderSpecific(Udot,VEC_CLASSID,4);
19033f72c5dSHong Zhang 
19133f72c5dSHong Zhang   ierr = PetscLogEventBegin(TS_JacobianPEval,ts,U,Amat,0);CHKERRQ(ierr);
19233f72c5dSHong Zhang   if (ts->ijacobianp) {
19333f72c5dSHong Zhang     PetscStackPush("TS user JacobianP function for sensitivity analysis");
19433f72c5dSHong Zhang     ierr = (*ts->ijacobianp)(ts,t,U,Udot,shift,Amat,ts->ijacobianpctx);CHKERRQ(ierr);
19533f72c5dSHong Zhang     PetscStackPop;
19633f72c5dSHong Zhang   }
19733f72c5dSHong Zhang   if (imex) {
19833f72c5dSHong Zhang     if (!ts->ijacobianp) {  /* system was written as Udot = G(t,U) */
19933f72c5dSHong Zhang       PetscBool assembled;
20033f72c5dSHong Zhang       ierr = MatZeroEntries(Amat);CHKERRQ(ierr);
20133f72c5dSHong Zhang       ierr = MatAssembled(Amat,&assembled);CHKERRQ(ierr);
20233f72c5dSHong Zhang       if (!assembled) {
20333f72c5dSHong Zhang         ierr = MatAssemblyBegin(Amat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
20433f72c5dSHong Zhang         ierr = MatAssemblyEnd(Amat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
20533f72c5dSHong Zhang       }
20633f72c5dSHong Zhang     }
20733f72c5dSHong Zhang   } else {
20833f72c5dSHong Zhang     if (ts->rhsjacobianp) {
20933f72c5dSHong Zhang       ierr = TSComputeRHSJacobianP(ts,t,U,ts->Jacprhs);CHKERRQ(ierr);
21033f72c5dSHong Zhang     }
21133f72c5dSHong Zhang     if (ts->Jacprhs == Amat) { /* No IJacobian, so we only have the RHS matrix */
21233f72c5dSHong Zhang       ierr = MatScale(Amat,-1);CHKERRQ(ierr);
21333f72c5dSHong Zhang     } else if (ts->Jacprhs) { /* Both IJacobian and RHSJacobian */
21433f72c5dSHong Zhang       MatStructure axpy = DIFFERENT_NONZERO_PATTERN;
21533f72c5dSHong Zhang       if (!ts->ijacobianp) { /* No IJacobianp provided, but we have a separate RHS matrix */
21633f72c5dSHong Zhang         ierr = MatZeroEntries(Amat);CHKERRQ(ierr);
21733f72c5dSHong Zhang       }
21833f72c5dSHong Zhang       ierr = MatAXPY(Amat,-1,ts->Jacprhs,axpy);CHKERRQ(ierr);
21933f72c5dSHong Zhang     }
22033f72c5dSHong Zhang   }
22133f72c5dSHong Zhang   ierr = PetscLogEventEnd(TS_JacobianPEval,ts,U,Amat,0);CHKERRQ(ierr);
22233f72c5dSHong Zhang   PetscFunctionReturn(0);
22333f72c5dSHong Zhang }
22433f72c5dSHong Zhang 
22533f72c5dSHong Zhang /*@C
226814e85d6SHong Zhang     TSSetCostIntegrand - Sets the routine for evaluating the integral term in one or more cost functions
227814e85d6SHong Zhang 
228814e85d6SHong Zhang     Logically Collective on TS
229814e85d6SHong Zhang 
230814e85d6SHong Zhang     Input Parameters:
231814e85d6SHong Zhang +   ts - the TS context obtained from TSCreate()
232814e85d6SHong Zhang .   numcost - number of gradients to be computed, this is the number of cost functions
233814e85d6SHong Zhang .   costintegral - vector that stores the integral values
234814e85d6SHong Zhang .   rf - routine for evaluating the integrand function
235c9ad9de0SHong Zhang .   drduf - function that computes the gradients of the r's with respect to u
236814e85d6SHong 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)
237814e85d6SHong Zhang .   fwd - flag indicating whether to evaluate cost integral in the forward run or the adjoint run
238814e85d6SHong Zhang -   ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
239814e85d6SHong Zhang 
240814e85d6SHong Zhang     Calling sequence of rf:
241c9ad9de0SHong Zhang $   PetscErrorCode rf(TS ts,PetscReal t,Vec U,Vec F,void *ctx);
242814e85d6SHong Zhang 
243c9ad9de0SHong Zhang     Calling sequence of drduf:
244c9ad9de0SHong Zhang $   PetscErroCode drduf(TS ts,PetscReal t,Vec U,Vec *dRdU,void *ctx);
245814e85d6SHong Zhang 
246814e85d6SHong Zhang     Calling sequence of drdpf:
247c9ad9de0SHong Zhang $   PetscErroCode drdpf(TS ts,PetscReal t,Vec U,Vec *dRdP,void *ctx);
248814e85d6SHong Zhang 
249cd4cee2dSHong Zhang     Level: deprecated
250814e85d6SHong Zhang 
251814e85d6SHong Zhang     Notes:
252814e85d6SHong Zhang     For optimization there is usually a single cost function (numcost = 1). For sensitivities there may be multiple cost functions
253814e85d6SHong Zhang 
254814e85d6SHong Zhang .seealso: TSSetRHSJacobianP(), TSGetCostGradients(), TSSetCostGradients()
255814e85d6SHong Zhang @*/
256814e85d6SHong Zhang PetscErrorCode TSSetCostIntegrand(TS ts,PetscInt numcost,Vec costintegral,PetscErrorCode (*rf)(TS,PetscReal,Vec,Vec,void*),
257c9ad9de0SHong Zhang                                                           PetscErrorCode (*drduf)(TS,PetscReal,Vec,Vec*,void*),
258814e85d6SHong Zhang                                                           PetscErrorCode (*drdpf)(TS,PetscReal,Vec,Vec*,void*),
259814e85d6SHong Zhang                                                           PetscBool fwd,void *ctx)
260814e85d6SHong Zhang {
261814e85d6SHong Zhang   PetscErrorCode ierr;
262814e85d6SHong Zhang 
263814e85d6SHong Zhang   PetscFunctionBegin;
264814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
265814e85d6SHong Zhang   if (costintegral) PetscValidHeaderSpecific(costintegral,VEC_CLASSID,3);
266814e85d6SHong 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()");
267814e85d6SHong Zhang   if (!ts->numcost) ts->numcost=numcost;
268814e85d6SHong Zhang 
269814e85d6SHong Zhang   if (costintegral) {
270814e85d6SHong Zhang     ierr = PetscObjectReference((PetscObject)costintegral);CHKERRQ(ierr);
271814e85d6SHong Zhang     ierr = VecDestroy(&ts->vec_costintegral);CHKERRQ(ierr);
272814e85d6SHong Zhang     ts->vec_costintegral = costintegral;
273814e85d6SHong Zhang   } else {
274814e85d6SHong Zhang     if (!ts->vec_costintegral) { /* Create a seq vec if user does not provide one */
275814e85d6SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,numcost,&ts->vec_costintegral);CHKERRQ(ierr);
276814e85d6SHong Zhang     } else {
277814e85d6SHong Zhang       ierr = VecSet(ts->vec_costintegral,0.0);CHKERRQ(ierr);
278814e85d6SHong Zhang     }
279814e85d6SHong Zhang   }
280814e85d6SHong Zhang   if (!ts->vec_costintegrand) {
281814e85d6SHong Zhang     ierr = VecDuplicate(ts->vec_costintegral,&ts->vec_costintegrand);CHKERRQ(ierr);
282814e85d6SHong Zhang   } else {
283814e85d6SHong Zhang     ierr = VecSet(ts->vec_costintegrand,0.0);CHKERRQ(ierr);
284814e85d6SHong Zhang   }
285814e85d6SHong Zhang   ts->costintegralfwd  = fwd; /* Evaluate the cost integral in forward run if fwd is true */
286814e85d6SHong Zhang   ts->costintegrand    = rf;
287814e85d6SHong Zhang   ts->costintegrandctx = ctx;
288c9ad9de0SHong Zhang   ts->drdufunction     = drduf;
289814e85d6SHong Zhang   ts->drdpfunction     = drdpf;
290814e85d6SHong Zhang   PetscFunctionReturn(0);
291814e85d6SHong Zhang }
292814e85d6SHong Zhang 
293b5b4017aSHong Zhang /*@C
294814e85d6SHong Zhang    TSGetCostIntegral - Returns the values of the integral term in the cost functions.
295814e85d6SHong Zhang    It is valid to call the routine after a backward run.
296814e85d6SHong Zhang 
297814e85d6SHong Zhang    Not Collective
298814e85d6SHong Zhang 
299814e85d6SHong Zhang    Input Parameter:
300814e85d6SHong Zhang .  ts - the TS context obtained from TSCreate()
301814e85d6SHong Zhang 
302814e85d6SHong Zhang    Output Parameter:
303814e85d6SHong Zhang .  v - the vector containing the integrals for each cost function
304814e85d6SHong Zhang 
305814e85d6SHong Zhang    Level: intermediate
306814e85d6SHong Zhang 
307814e85d6SHong Zhang .seealso: TSSetCostIntegrand()
308814e85d6SHong Zhang 
309814e85d6SHong Zhang @*/
310814e85d6SHong Zhang PetscErrorCode  TSGetCostIntegral(TS ts,Vec *v)
311814e85d6SHong Zhang {
312cd4cee2dSHong Zhang   TS             quadts;
313cd4cee2dSHong Zhang   PetscErrorCode ierr;
314cd4cee2dSHong Zhang 
315814e85d6SHong Zhang   PetscFunctionBegin;
316814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
317814e85d6SHong Zhang   PetscValidPointer(v,2);
318cd4cee2dSHong Zhang   ierr = TSGetQuadratureTS(ts,NULL,&quadts);CHKERRQ(ierr);
319cd4cee2dSHong Zhang   *v = quadts->vec_sol;
320814e85d6SHong Zhang   PetscFunctionReturn(0);
321814e85d6SHong Zhang }
322814e85d6SHong Zhang 
323b5b4017aSHong Zhang /*@C
324814e85d6SHong Zhang    TSComputeCostIntegrand - Evaluates the integral function in the cost functions.
325814e85d6SHong Zhang 
326814e85d6SHong Zhang    Input Parameters:
327814e85d6SHong Zhang +  ts - the TS context
328814e85d6SHong Zhang .  t - current time
329c9ad9de0SHong Zhang -  U - state vector, i.e. current solution
330814e85d6SHong Zhang 
331814e85d6SHong Zhang    Output Parameter:
332c9ad9de0SHong Zhang .  Q - vector of size numcost to hold the outputs
333814e85d6SHong Zhang 
334369cf35fSHong Zhang    Notes:
335814e85d6SHong Zhang    Most users should not need to explicitly call this routine, as it
336814e85d6SHong Zhang    is used internally within the sensitivity analysis context.
337814e85d6SHong Zhang 
338cd4cee2dSHong Zhang    Level: deprecated
339814e85d6SHong Zhang 
340814e85d6SHong Zhang .seealso: TSSetCostIntegrand()
341814e85d6SHong Zhang @*/
342c9ad9de0SHong Zhang PetscErrorCode TSComputeCostIntegrand(TS ts,PetscReal t,Vec U,Vec Q)
343814e85d6SHong Zhang {
344814e85d6SHong Zhang   PetscErrorCode ierr;
345814e85d6SHong Zhang 
346814e85d6SHong Zhang   PetscFunctionBegin;
347814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
348c9ad9de0SHong Zhang   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
349c9ad9de0SHong Zhang   PetscValidHeaderSpecific(Q,VEC_CLASSID,4);
350814e85d6SHong Zhang 
351c9ad9de0SHong Zhang   ierr = PetscLogEventBegin(TS_FunctionEval,ts,U,Q,0);CHKERRQ(ierr);
352814e85d6SHong Zhang   if (ts->costintegrand) {
353814e85d6SHong Zhang     PetscStackPush("TS user integrand in the cost function");
354c9ad9de0SHong Zhang     ierr = (*ts->costintegrand)(ts,t,U,Q,ts->costintegrandctx);CHKERRQ(ierr);
355814e85d6SHong Zhang     PetscStackPop;
356814e85d6SHong Zhang   } else {
357c9ad9de0SHong Zhang     ierr = VecZeroEntries(Q);CHKERRQ(ierr);
358814e85d6SHong Zhang   }
359814e85d6SHong Zhang 
360c9ad9de0SHong Zhang   ierr = PetscLogEventEnd(TS_FunctionEval,ts,U,Q,0);CHKERRQ(ierr);
361814e85d6SHong Zhang   PetscFunctionReturn(0);
362814e85d6SHong Zhang }
363814e85d6SHong Zhang 
364b5b4017aSHong Zhang /*@C
365d76be551SHong Zhang   TSComputeDRDUFunction - Deprecated, use TSGetQuadratureTS() then TSComputeRHSJacobian()
366814e85d6SHong Zhang 
367d76be551SHong Zhang   Level: deprecated
368814e85d6SHong Zhang 
369814e85d6SHong Zhang @*/
370c9ad9de0SHong Zhang PetscErrorCode TSComputeDRDUFunction(TS ts,PetscReal t,Vec U,Vec *DRDU)
371814e85d6SHong Zhang {
372814e85d6SHong Zhang   PetscErrorCode ierr;
373814e85d6SHong Zhang 
374814e85d6SHong Zhang   PetscFunctionBegin;
37533f72c5dSHong Zhang   if (!DRDU) PetscFunctionReturn(0);
376814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
377c9ad9de0SHong Zhang   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
378814e85d6SHong Zhang 
379c9ad9de0SHong Zhang   PetscStackPush("TS user DRDU function for sensitivity analysis");
380c9ad9de0SHong Zhang   ierr = (*ts->drdufunction)(ts,t,U,DRDU,ts->costintegrandctx); CHKERRQ(ierr);
381814e85d6SHong Zhang   PetscStackPop;
382814e85d6SHong Zhang   PetscFunctionReturn(0);
383814e85d6SHong Zhang }
384814e85d6SHong Zhang 
385b5b4017aSHong Zhang /*@C
386d76be551SHong Zhang   TSComputeDRDPFunction - Deprecated, use TSGetQuadratureTS() then TSComputeRHSJacobianP()
387814e85d6SHong Zhang 
388d76be551SHong Zhang   Level: deprecated
389814e85d6SHong Zhang 
390814e85d6SHong Zhang @*/
391c9ad9de0SHong Zhang PetscErrorCode TSComputeDRDPFunction(TS ts,PetscReal t,Vec U,Vec *DRDP)
392814e85d6SHong Zhang {
393814e85d6SHong Zhang   PetscErrorCode ierr;
394814e85d6SHong Zhang 
395814e85d6SHong Zhang   PetscFunctionBegin;
39633f72c5dSHong Zhang   if (!DRDP) PetscFunctionReturn(0);
397814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
398c9ad9de0SHong Zhang   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
399814e85d6SHong Zhang 
400814e85d6SHong Zhang   PetscStackPush("TS user DRDP function for sensitivity analysis");
401c9ad9de0SHong Zhang   ierr = (*ts->drdpfunction)(ts,t,U,DRDP,ts->costintegrandctx); CHKERRQ(ierr);
402814e85d6SHong Zhang   PetscStackPop;
403814e85d6SHong Zhang   PetscFunctionReturn(0);
404814e85d6SHong Zhang }
405814e85d6SHong Zhang 
406b5b4017aSHong Zhang /*@C
407b5b4017aSHong 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.
408b5b4017aSHong Zhang 
409b5b4017aSHong Zhang   Logically Collective on TS
410b5b4017aSHong Zhang 
411b5b4017aSHong Zhang   Input Parameters:
412b5b4017aSHong Zhang + ts - TS context obtained from TSCreate()
413b5b4017aSHong Zhang . ihp1 - an array of vectors storing the result of vector-Hessian-vector product for F_UU
414b5b4017aSHong Zhang . hessianproductfunc1 - vector-Hessian-vector product function for F_UU
415b5b4017aSHong Zhang . ihp2 - an array of vectors storing the result of vector-Hessian-vector product for F_UP
416b5b4017aSHong Zhang . hessianproductfunc2 - vector-Hessian-vector product function for F_UP
417b5b4017aSHong Zhang . ihp3 - an array of vectors storing the result of vector-Hessian-vector product for F_PU
418b5b4017aSHong Zhang . hessianproductfunc3 - vector-Hessian-vector product function for F_PU
419b5b4017aSHong Zhang . ihp4 - an array of vectors storing the result of vector-Hessian-vector product for F_PP
420b5b4017aSHong Zhang . hessianproductfunc4 - vector-Hessian-vector product function for F_PP
421b5b4017aSHong Zhang 
422b5b4017aSHong Zhang   Calling sequence of ihessianproductfunc:
423369cf35fSHong Zhang $ ihessianproductfunc (TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx);
424b5b4017aSHong Zhang +   t - current timestep
425b5b4017aSHong Zhang .   U - input vector (current ODE solution)
426369cf35fSHong Zhang .   Vl - an array of input vectors to be left-multiplied with the Hessian
427b5b4017aSHong Zhang .   Vr - input vector to be right-multiplied with the Hessian
428369cf35fSHong Zhang .   VHV - an array of output vectors for vector-Hessian-vector product
429b5b4017aSHong Zhang -   ctx - [optional] user-defined function context
430b5b4017aSHong Zhang 
431b5b4017aSHong Zhang   Level: intermediate
432b5b4017aSHong Zhang 
433369cf35fSHong Zhang   Notes:
434369cf35fSHong Zhang   The first Hessian function and the working array are required.
435369cf35fSHong Zhang   As an example to implement the callback functions, the second callback function calculates the vector-Hessian-vector product
436369cf35fSHong Zhang   $ Vl_n^T*F_UP*Vr
437369cf35fSHong 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.
438369cf35fSHong Zhang   Each entry of F_UP corresponds to the derivative
439369cf35fSHong Zhang   $ F_UP[i][j][k] = \frac{\partial^2 F[i]}{\partial U[j] \partial P[k]}.
440369cf35fSHong 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
441369cf35fSHong Zhang   $ VHV_n[j] = \sum_i \sum_k {Vl_n[i] * F_UP[i][j][k] * Vr[k]}
442369cf35fSHong Zhang   If the cost function is a scalar, there will be only one vector in Vl and VHV.
443b5b4017aSHong Zhang 
444b5b4017aSHong Zhang .seealso:
445b5b4017aSHong Zhang @*/
446b5b4017aSHong Zhang PetscErrorCode TSSetIHessianProduct(TS ts,Vec *ihp1,PetscErrorCode (*ihessianproductfunc1)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*),
447b5b4017aSHong Zhang                                           Vec *ihp2,PetscErrorCode (*ihessianproductfunc2)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*),
448b5b4017aSHong Zhang                                           Vec *ihp3,PetscErrorCode (*ihessianproductfunc3)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*),
449b5b4017aSHong Zhang                                           Vec *ihp4,PetscErrorCode (*ihessianproductfunc4)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*),
450b5b4017aSHong Zhang                                     void *ctx)
451b5b4017aSHong Zhang {
452b5b4017aSHong Zhang   PetscFunctionBegin;
453b5b4017aSHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
454b5b4017aSHong Zhang   PetscValidPointer(ihp1,2);
455b5b4017aSHong Zhang 
456b5b4017aSHong Zhang   ts->ihessianproductctx = ctx;
457b5b4017aSHong Zhang   if (ihp1) ts->vecs_fuu = ihp1;
458b5b4017aSHong Zhang   if (ihp2) ts->vecs_fup = ihp2;
459b5b4017aSHong Zhang   if (ihp3) ts->vecs_fpu = ihp3;
460b5b4017aSHong Zhang   if (ihp4) ts->vecs_fpp = ihp4;
461b5b4017aSHong Zhang   ts->ihessianproduct_fuu = ihessianproductfunc1;
462b5b4017aSHong Zhang   ts->ihessianproduct_fup = ihessianproductfunc2;
463b5b4017aSHong Zhang   ts->ihessianproduct_fpu = ihessianproductfunc3;
464b5b4017aSHong Zhang   ts->ihessianproduct_fpp = ihessianproductfunc4;
465b5b4017aSHong Zhang   PetscFunctionReturn(0);
466b5b4017aSHong Zhang }
467b5b4017aSHong Zhang 
468b5b4017aSHong Zhang /*@C
469b1e111ebSHong Zhang   TSComputeIHessianProductFunctionUU - Runs the user-defined vector-Hessian-vector product function for Fuu.
470b5b4017aSHong Zhang 
471b5b4017aSHong Zhang   Collective on TS
472b5b4017aSHong Zhang 
473b5b4017aSHong Zhang   Input Parameters:
474b5b4017aSHong Zhang . ts   - The TS context obtained from TSCreate()
475b5b4017aSHong Zhang 
476b5b4017aSHong Zhang   Notes:
477b1e111ebSHong Zhang   TSComputeIHessianProductFunctionUU() is typically used for sensitivity implementation,
478b5b4017aSHong Zhang   so most users would not generally call this routine themselves.
479b5b4017aSHong Zhang 
480b5b4017aSHong Zhang   Level: developer
481b5b4017aSHong Zhang 
482b5b4017aSHong Zhang .seealso: TSSetIHessianProduct()
483b5b4017aSHong Zhang @*/
484b1e111ebSHong Zhang PetscErrorCode TSComputeIHessianProductFunctionUU(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV)
485b5b4017aSHong Zhang {
486b5b4017aSHong Zhang   PetscErrorCode ierr;
487b5b4017aSHong Zhang 
488b5b4017aSHong Zhang   PetscFunctionBegin;
48933f72c5dSHong Zhang   if (!VHV) PetscFunctionReturn(0);
490b5b4017aSHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
491b5b4017aSHong Zhang   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
492b5b4017aSHong Zhang 
49367633408SHong Zhang   if (ts->ihessianproduct_fuu) {
494b5b4017aSHong Zhang     PetscStackPush("TS user IHessianProduct function 1 for sensitivity analysis");
495b5b4017aSHong Zhang     ierr = (*ts->ihessianproduct_fuu)(ts,t,U,Vl,Vr,VHV,ts->ihessianproductctx);CHKERRQ(ierr);
496b5b4017aSHong Zhang     PetscStackPop;
49767633408SHong Zhang   }
49867633408SHong Zhang   /* does not consider IMEX for now, so either IHessian or RHSHessian will be calculated, using the same output VHV */
49967633408SHong Zhang   if (ts->rhshessianproduct_guu) {
50067633408SHong Zhang     PetscInt nadj;
501b1e111ebSHong Zhang     ierr = TSComputeRHSHessianProductFunctionUU(ts,t,U,Vl,Vr,VHV);CHKERRQ(ierr);
50267633408SHong Zhang     for (nadj=0; nadj<ts->numcost; nadj++) {
50367633408SHong Zhang       ierr = VecScale(VHV[nadj],-1);CHKERRQ(ierr);
50467633408SHong Zhang     }
50567633408SHong Zhang   }
506b5b4017aSHong Zhang   PetscFunctionReturn(0);
507b5b4017aSHong Zhang }
508b5b4017aSHong Zhang 
509b5b4017aSHong Zhang /*@C
510b1e111ebSHong Zhang   TSComputeIHessianProductFunctionUP - Runs the user-defined vector-Hessian-vector product function for Fup.
511b5b4017aSHong Zhang 
512b5b4017aSHong Zhang   Collective on TS
513b5b4017aSHong Zhang 
514b5b4017aSHong Zhang   Input Parameters:
515b5b4017aSHong Zhang . ts   - The TS context obtained from TSCreate()
516b5b4017aSHong Zhang 
517b5b4017aSHong Zhang   Notes:
518b1e111ebSHong Zhang   TSComputeIHessianProductFunctionUP() is typically used for sensitivity implementation,
519b5b4017aSHong Zhang   so most users would not generally call this routine themselves.
520b5b4017aSHong Zhang 
521b5b4017aSHong Zhang   Level: developer
522b5b4017aSHong Zhang 
523b5b4017aSHong Zhang .seealso: TSSetIHessianProduct()
524b5b4017aSHong Zhang @*/
525b1e111ebSHong Zhang PetscErrorCode TSComputeIHessianProductFunctionUP(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV)
526b5b4017aSHong Zhang {
527b5b4017aSHong Zhang   PetscErrorCode ierr;
528b5b4017aSHong Zhang 
529b5b4017aSHong Zhang   PetscFunctionBegin;
53033f72c5dSHong Zhang   if (!VHV) PetscFunctionReturn(0);
531b5b4017aSHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
532b5b4017aSHong Zhang   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
533b5b4017aSHong Zhang 
53467633408SHong Zhang   if (ts->ihessianproduct_fup) {
535b5b4017aSHong Zhang     PetscStackPush("TS user IHessianProduct function 2 for sensitivity analysis");
536b5b4017aSHong Zhang     ierr = (*ts->ihessianproduct_fup)(ts,t,U,Vl,Vr,VHV,ts->ihessianproductctx);CHKERRQ(ierr);
537b5b4017aSHong Zhang     PetscStackPop;
53867633408SHong Zhang   }
53967633408SHong Zhang   /* does not consider IMEX for now, so either IHessian or RHSHessian will be calculated, using the same output VHV */
54067633408SHong Zhang   if (ts->rhshessianproduct_gup) {
54167633408SHong Zhang     PetscInt nadj;
542b1e111ebSHong Zhang     ierr = TSComputeRHSHessianProductFunctionUP(ts,t,U,Vl,Vr,VHV);CHKERRQ(ierr);
54367633408SHong Zhang     for (nadj=0; nadj<ts->numcost; nadj++) {
54467633408SHong Zhang       ierr = VecScale(VHV[nadj],-1);CHKERRQ(ierr);
54567633408SHong Zhang     }
54667633408SHong Zhang   }
547b5b4017aSHong Zhang   PetscFunctionReturn(0);
548b5b4017aSHong Zhang }
549b5b4017aSHong Zhang 
550b5b4017aSHong Zhang /*@C
551b1e111ebSHong Zhang   TSComputeIHessianProductFunctionPU - Runs the user-defined vector-Hessian-vector product function for Fpu.
552b5b4017aSHong Zhang 
553b5b4017aSHong Zhang   Collective on TS
554b5b4017aSHong Zhang 
555b5b4017aSHong Zhang   Input Parameters:
556b5b4017aSHong Zhang . ts   - The TS context obtained from TSCreate()
557b5b4017aSHong Zhang 
558b5b4017aSHong Zhang   Notes:
559b1e111ebSHong Zhang   TSComputeIHessianProductFunctionPU() is typically used for sensitivity implementation,
560b5b4017aSHong Zhang   so most users would not generally call this routine themselves.
561b5b4017aSHong Zhang 
562b5b4017aSHong Zhang   Level: developer
563b5b4017aSHong Zhang 
564b5b4017aSHong Zhang .seealso: TSSetIHessianProduct()
565b5b4017aSHong Zhang @*/
566b1e111ebSHong Zhang PetscErrorCode TSComputeIHessianProductFunctionPU(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV)
567b5b4017aSHong Zhang {
568b5b4017aSHong Zhang   PetscErrorCode ierr;
569b5b4017aSHong Zhang 
570b5b4017aSHong Zhang   PetscFunctionBegin;
57133f72c5dSHong Zhang   if (!VHV) PetscFunctionReturn(0);
572b5b4017aSHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
573b5b4017aSHong Zhang   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
574b5b4017aSHong Zhang 
57567633408SHong Zhang   if (ts->ihessianproduct_fpu) {
576b5b4017aSHong Zhang     PetscStackPush("TS user IHessianProduct function 3 for sensitivity analysis");
577b5b4017aSHong Zhang     ierr = (*ts->ihessianproduct_fpu)(ts,t,U,Vl,Vr,VHV,ts->ihessianproductctx);CHKERRQ(ierr);
578b5b4017aSHong Zhang     PetscStackPop;
57967633408SHong Zhang   }
58067633408SHong Zhang   /* does not consider IMEX for now, so either IHessian or RHSHessian will be calculated, using the same output VHV */
58167633408SHong Zhang   if (ts->rhshessianproduct_gpu) {
58267633408SHong Zhang     PetscInt nadj;
583b1e111ebSHong Zhang     ierr = TSComputeRHSHessianProductFunctionPU(ts,t,U,Vl,Vr,VHV);CHKERRQ(ierr);
58467633408SHong Zhang     for (nadj=0; nadj<ts->numcost; nadj++) {
58567633408SHong Zhang       ierr = VecScale(VHV[nadj],-1);CHKERRQ(ierr);
58667633408SHong Zhang     }
58767633408SHong Zhang   }
588b5b4017aSHong Zhang   PetscFunctionReturn(0);
589b5b4017aSHong Zhang }
590b5b4017aSHong Zhang 
591b5b4017aSHong Zhang /*@C
592b1e111ebSHong Zhang   TSComputeIHessianProductFunctionPP - Runs the user-defined vector-Hessian-vector product function for Fpp.
593b5b4017aSHong Zhang 
594b5b4017aSHong Zhang   Collective on TS
595b5b4017aSHong Zhang 
596b5b4017aSHong Zhang   Input Parameters:
597b5b4017aSHong Zhang . ts   - The TS context obtained from TSCreate()
598b5b4017aSHong Zhang 
599b5b4017aSHong Zhang   Notes:
600b1e111ebSHong Zhang   TSComputeIHessianProductFunctionPP() is typically used for sensitivity implementation,
601b5b4017aSHong Zhang   so most users would not generally call this routine themselves.
602b5b4017aSHong Zhang 
603b5b4017aSHong Zhang   Level: developer
604b5b4017aSHong Zhang 
605b5b4017aSHong Zhang .seealso: TSSetIHessianProduct()
606b5b4017aSHong Zhang @*/
607b1e111ebSHong Zhang PetscErrorCode TSComputeIHessianProductFunctionPP(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV)
608b5b4017aSHong Zhang {
609b5b4017aSHong Zhang   PetscErrorCode ierr;
610b5b4017aSHong Zhang 
611b5b4017aSHong Zhang   PetscFunctionBegin;
61233f72c5dSHong Zhang   if (!VHV) PetscFunctionReturn(0);
613b5b4017aSHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
614b5b4017aSHong Zhang   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
615b5b4017aSHong Zhang 
61667633408SHong Zhang   if (ts->ihessianproduct_fpp) {
617b5b4017aSHong Zhang     PetscStackPush("TS user IHessianProduct function 3 for sensitivity analysis");
618b5b4017aSHong Zhang     ierr = (*ts->ihessianproduct_fpp)(ts,t,U,Vl,Vr,VHV,ts->ihessianproductctx);CHKERRQ(ierr);
619b5b4017aSHong Zhang     PetscStackPop;
62067633408SHong Zhang   }
62167633408SHong Zhang   /* does not consider IMEX for now, so either IHessian or RHSHessian will be calculated, using the same output VHV */
62267633408SHong Zhang   if (ts->rhshessianproduct_gpp) {
62367633408SHong Zhang     PetscInt nadj;
624b1e111ebSHong Zhang     ierr = TSComputeRHSHessianProductFunctionPP(ts,t,U,Vl,Vr,VHV);CHKERRQ(ierr);
62567633408SHong Zhang     for (nadj=0; nadj<ts->numcost; nadj++) {
62667633408SHong Zhang       ierr = VecScale(VHV[nadj],-1);CHKERRQ(ierr);
62767633408SHong Zhang     }
62867633408SHong Zhang   }
629b5b4017aSHong Zhang   PetscFunctionReturn(0);
630b5b4017aSHong Zhang }
631b5b4017aSHong Zhang 
63213af1a74SHong Zhang /*@C
6334c500f23SPierre Jolivet   TSSetRHSHessianProduct - Sets the function that computes the vector-Hessian-vector product. The Hessian is the second-order derivative of G (RHSFunction) w.r.t. the state variable.
63413af1a74SHong Zhang 
63513af1a74SHong Zhang   Logically Collective on TS
63613af1a74SHong Zhang 
63713af1a74SHong Zhang   Input Parameters:
63813af1a74SHong Zhang + ts - TS context obtained from TSCreate()
63913af1a74SHong Zhang . rhshp1 - an array of vectors storing the result of vector-Hessian-vector product for G_UU
64013af1a74SHong Zhang . hessianproductfunc1 - vector-Hessian-vector product function for G_UU
64113af1a74SHong Zhang . rhshp2 - an array of vectors storing the result of vector-Hessian-vector product for G_UP
64213af1a74SHong Zhang . hessianproductfunc2 - vector-Hessian-vector product function for G_UP
64313af1a74SHong Zhang . rhshp3 - an array of vectors storing the result of vector-Hessian-vector product for G_PU
64413af1a74SHong Zhang . hessianproductfunc3 - vector-Hessian-vector product function for G_PU
64513af1a74SHong Zhang . rhshp4 - an array of vectors storing the result of vector-Hessian-vector product for G_PP
64613af1a74SHong Zhang . hessianproductfunc4 - vector-Hessian-vector product function for G_PP
64713af1a74SHong Zhang 
64813af1a74SHong Zhang   Calling sequence of ihessianproductfunc:
649369cf35fSHong Zhang $ rhshessianproductfunc (TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx);
65013af1a74SHong Zhang +   t - current timestep
65113af1a74SHong Zhang .   U - input vector (current ODE solution)
652369cf35fSHong Zhang .   Vl - an array of input vectors to be left-multiplied with the Hessian
65313af1a74SHong Zhang .   Vr - input vector to be right-multiplied with the Hessian
654369cf35fSHong Zhang .   VHV - an array of output vectors for vector-Hessian-vector product
65513af1a74SHong Zhang -   ctx - [optional] user-defined function context
65613af1a74SHong Zhang 
65713af1a74SHong Zhang   Level: intermediate
65813af1a74SHong Zhang 
659369cf35fSHong Zhang   Notes:
660369cf35fSHong Zhang   The first Hessian function and the working array are required.
661369cf35fSHong Zhang   As an example to implement the callback functions, the second callback function calculates the vector-Hessian-vector product
662369cf35fSHong Zhang   $ Vl_n^T*G_UP*Vr
663369cf35fSHong 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.
664369cf35fSHong Zhang   Each entry of G_UP corresponds to the derivative
665369cf35fSHong Zhang   $ G_UP[i][j][k] = \frac{\partial^2 G[i]}{\partial U[j] \partial P[k]}.
666369cf35fSHong 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
667369cf35fSHong Zhang   $ VHV_n[j] = \sum_i \sum_k {Vl_n[i] * G_UP[i][j][k] * Vr[k]}
668369cf35fSHong Zhang   If the cost function is a scalar, there will be only one vector in Vl and VHV.
66913af1a74SHong Zhang 
67013af1a74SHong Zhang .seealso:
67113af1a74SHong Zhang @*/
67213af1a74SHong Zhang PetscErrorCode TSSetRHSHessianProduct(TS ts,Vec *rhshp1,PetscErrorCode (*rhshessianproductfunc1)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*),
67313af1a74SHong Zhang                                           Vec *rhshp2,PetscErrorCode (*rhshessianproductfunc2)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*),
67413af1a74SHong Zhang                                           Vec *rhshp3,PetscErrorCode (*rhshessianproductfunc3)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*),
67513af1a74SHong Zhang                                           Vec *rhshp4,PetscErrorCode (*rhshessianproductfunc4)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*),
67613af1a74SHong Zhang                                     void *ctx)
67713af1a74SHong Zhang {
67813af1a74SHong Zhang   PetscFunctionBegin;
67913af1a74SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
68013af1a74SHong Zhang   PetscValidPointer(rhshp1,2);
68113af1a74SHong Zhang 
68213af1a74SHong Zhang   ts->rhshessianproductctx = ctx;
68313af1a74SHong Zhang   if (rhshp1) ts->vecs_guu = rhshp1;
68413af1a74SHong Zhang   if (rhshp2) ts->vecs_gup = rhshp2;
68513af1a74SHong Zhang   if (rhshp3) ts->vecs_gpu = rhshp3;
68613af1a74SHong Zhang   if (rhshp4) ts->vecs_gpp = rhshp4;
68713af1a74SHong Zhang   ts->rhshessianproduct_guu = rhshessianproductfunc1;
68813af1a74SHong Zhang   ts->rhshessianproduct_gup = rhshessianproductfunc2;
68913af1a74SHong Zhang   ts->rhshessianproduct_gpu = rhshessianproductfunc3;
69013af1a74SHong Zhang   ts->rhshessianproduct_gpp = rhshessianproductfunc4;
69113af1a74SHong Zhang   PetscFunctionReturn(0);
69213af1a74SHong Zhang }
69313af1a74SHong Zhang 
69413af1a74SHong Zhang /*@C
695b1e111ebSHong Zhang   TSComputeRHSHessianProductFunctionUU - Runs the user-defined vector-Hessian-vector product function for Guu.
69613af1a74SHong Zhang 
69713af1a74SHong Zhang   Collective on TS
69813af1a74SHong Zhang 
69913af1a74SHong Zhang   Input Parameters:
70013af1a74SHong Zhang . ts   - The TS context obtained from TSCreate()
70113af1a74SHong Zhang 
70213af1a74SHong Zhang   Notes:
703b1e111ebSHong Zhang   TSComputeRHSHessianProductFunctionUU() is typically used for sensitivity implementation,
70413af1a74SHong Zhang   so most users would not generally call this routine themselves.
70513af1a74SHong Zhang 
70613af1a74SHong Zhang   Level: developer
70713af1a74SHong Zhang 
70813af1a74SHong Zhang .seealso: TSSetRHSHessianProduct()
70913af1a74SHong Zhang @*/
710b1e111ebSHong Zhang PetscErrorCode TSComputeRHSHessianProductFunctionUU(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV)
71113af1a74SHong Zhang {
71213af1a74SHong Zhang   PetscErrorCode ierr;
71313af1a74SHong Zhang 
71413af1a74SHong Zhang   PetscFunctionBegin;
71513af1a74SHong Zhang   if (!VHV) PetscFunctionReturn(0);
71613af1a74SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
71713af1a74SHong Zhang   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
71813af1a74SHong Zhang 
71913af1a74SHong Zhang   PetscStackPush("TS user RHSHessianProduct function 1 for sensitivity analysis");
72013af1a74SHong Zhang   ierr = (*ts->rhshessianproduct_guu)(ts,t,U,Vl,Vr,VHV,ts->rhshessianproductctx);CHKERRQ(ierr);
72113af1a74SHong Zhang   PetscStackPop;
72213af1a74SHong Zhang   PetscFunctionReturn(0);
72313af1a74SHong Zhang }
72413af1a74SHong Zhang 
72513af1a74SHong Zhang /*@C
726b1e111ebSHong Zhang   TSComputeRHSHessianProductFunctionUP - Runs the user-defined vector-Hessian-vector product function for Gup.
72713af1a74SHong Zhang 
72813af1a74SHong Zhang   Collective on TS
72913af1a74SHong Zhang 
73013af1a74SHong Zhang   Input Parameters:
73113af1a74SHong Zhang . ts   - The TS context obtained from TSCreate()
73213af1a74SHong Zhang 
73313af1a74SHong Zhang   Notes:
734b1e111ebSHong Zhang   TSComputeRHSHessianProductFunctionUP() is typically used for sensitivity implementation,
73513af1a74SHong Zhang   so most users would not generally call this routine themselves.
73613af1a74SHong Zhang 
73713af1a74SHong Zhang   Level: developer
73813af1a74SHong Zhang 
73913af1a74SHong Zhang .seealso: TSSetRHSHessianProduct()
74013af1a74SHong Zhang @*/
741b1e111ebSHong Zhang PetscErrorCode TSComputeRHSHessianProductFunctionUP(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV)
74213af1a74SHong Zhang {
74313af1a74SHong Zhang   PetscErrorCode ierr;
74413af1a74SHong Zhang 
74513af1a74SHong Zhang   PetscFunctionBegin;
74613af1a74SHong Zhang   if (!VHV) PetscFunctionReturn(0);
74713af1a74SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
74813af1a74SHong Zhang   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
74913af1a74SHong Zhang 
75013af1a74SHong Zhang   PetscStackPush("TS user RHSHessianProduct function 2 for sensitivity analysis");
75113af1a74SHong Zhang   ierr = (*ts->rhshessianproduct_gup)(ts,t,U,Vl,Vr,VHV,ts->rhshessianproductctx);CHKERRQ(ierr);
75213af1a74SHong Zhang   PetscStackPop;
75313af1a74SHong Zhang   PetscFunctionReturn(0);
75413af1a74SHong Zhang }
75513af1a74SHong Zhang 
75613af1a74SHong Zhang /*@C
757b1e111ebSHong Zhang   TSComputeRHSHessianProductFunctionPU - Runs the user-defined vector-Hessian-vector product function for Gpu.
75813af1a74SHong Zhang 
75913af1a74SHong Zhang   Collective on TS
76013af1a74SHong Zhang 
76113af1a74SHong Zhang   Input Parameters:
76213af1a74SHong Zhang . ts   - The TS context obtained from TSCreate()
76313af1a74SHong Zhang 
76413af1a74SHong Zhang   Notes:
765b1e111ebSHong Zhang   TSComputeRHSHessianProductFunctionPU() is typically used for sensitivity implementation,
76613af1a74SHong Zhang   so most users would not generally call this routine themselves.
76713af1a74SHong Zhang 
76813af1a74SHong Zhang   Level: developer
76913af1a74SHong Zhang 
77013af1a74SHong Zhang .seealso: TSSetRHSHessianProduct()
77113af1a74SHong Zhang @*/
772b1e111ebSHong Zhang PetscErrorCode TSComputeRHSHessianProductFunctionPU(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV)
77313af1a74SHong Zhang {
77413af1a74SHong Zhang   PetscErrorCode ierr;
77513af1a74SHong Zhang 
77613af1a74SHong Zhang   PetscFunctionBegin;
77713af1a74SHong Zhang   if (!VHV) PetscFunctionReturn(0);
77813af1a74SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
77913af1a74SHong Zhang   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
78013af1a74SHong Zhang 
78113af1a74SHong Zhang   PetscStackPush("TS user RHSHessianProduct function 3 for sensitivity analysis");
78213af1a74SHong Zhang   ierr = (*ts->rhshessianproduct_gpu)(ts,t,U,Vl,Vr,VHV,ts->rhshessianproductctx);CHKERRQ(ierr);
78313af1a74SHong Zhang   PetscStackPop;
78413af1a74SHong Zhang   PetscFunctionReturn(0);
78513af1a74SHong Zhang }
78613af1a74SHong Zhang 
78713af1a74SHong Zhang /*@C
788b1e111ebSHong Zhang   TSComputeRHSHessianProductFunctionPP - Runs the user-defined vector-Hessian-vector product function for Gpp.
78913af1a74SHong Zhang 
79013af1a74SHong Zhang   Collective on TS
79113af1a74SHong Zhang 
79213af1a74SHong Zhang   Input Parameters:
79313af1a74SHong Zhang . ts   - The TS context obtained from TSCreate()
79413af1a74SHong Zhang 
79513af1a74SHong Zhang   Notes:
796b1e111ebSHong Zhang   TSComputeRHSHessianProductFunctionPP() is typically used for sensitivity implementation,
79713af1a74SHong Zhang   so most users would not generally call this routine themselves.
79813af1a74SHong Zhang 
79913af1a74SHong Zhang   Level: developer
80013af1a74SHong Zhang 
80113af1a74SHong Zhang .seealso: TSSetRHSHessianProduct()
80213af1a74SHong Zhang @*/
803b1e111ebSHong Zhang PetscErrorCode TSComputeRHSHessianProductFunctionPP(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV)
80413af1a74SHong Zhang {
80513af1a74SHong Zhang   PetscErrorCode ierr;
80613af1a74SHong Zhang 
80713af1a74SHong Zhang   PetscFunctionBegin;
80813af1a74SHong Zhang   if (!VHV) PetscFunctionReturn(0);
80913af1a74SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
81013af1a74SHong Zhang   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
81113af1a74SHong Zhang 
81213af1a74SHong Zhang   PetscStackPush("TS user RHSHessianProduct function 3 for sensitivity analysis");
81313af1a74SHong Zhang   ierr = (*ts->rhshessianproduct_gpp)(ts,t,U,Vl,Vr,VHV,ts->rhshessianproductctx);CHKERRQ(ierr);
81413af1a74SHong Zhang   PetscStackPop;
81513af1a74SHong Zhang   PetscFunctionReturn(0);
81613af1a74SHong Zhang }
81713af1a74SHong Zhang 
818814e85d6SHong Zhang /* --------------------------- Adjoint sensitivity ---------------------------*/
819814e85d6SHong Zhang 
820814e85d6SHong Zhang /*@
821814e85d6SHong 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
822814e85d6SHong Zhang       for use by the TSAdjoint routines.
823814e85d6SHong Zhang 
824d083f849SBarry Smith    Logically Collective on TS
825814e85d6SHong Zhang 
826814e85d6SHong Zhang    Input Parameters:
827814e85d6SHong Zhang +  ts - the TS context obtained from TSCreate()
828814e85d6SHong 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
829814e85d6SHong Zhang -  mu - gradients with respect to the parameters, the number of entries in these vectors is the same as the number of parameters
830814e85d6SHong Zhang 
831814e85d6SHong Zhang    Level: beginner
832814e85d6SHong Zhang 
833814e85d6SHong Zhang    Notes:
834814e85d6SHong Zhang     the entries in these vectors must be correctly initialized with the values lamda_i = df/dy|finaltime  mu_i = df/dp|finaltime
835814e85d6SHong Zhang 
836814e85d6SHong Zhang    After TSAdjointSolve() is called the lamba and the mu contain the computed sensitivities
837814e85d6SHong Zhang 
838b5b4017aSHong Zhang .seealso TSGetCostGradients()
839814e85d6SHong Zhang @*/
840814e85d6SHong Zhang PetscErrorCode TSSetCostGradients(TS ts,PetscInt numcost,Vec *lambda,Vec *mu)
841814e85d6SHong Zhang {
842814e85d6SHong Zhang   PetscFunctionBegin;
843814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
844814e85d6SHong Zhang   PetscValidPointer(lambda,2);
845814e85d6SHong Zhang   ts->vecs_sensi  = lambda;
846814e85d6SHong Zhang   ts->vecs_sensip = mu;
847814e85d6SHong 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");
848814e85d6SHong Zhang   ts->numcost  = numcost;
849814e85d6SHong Zhang   PetscFunctionReturn(0);
850814e85d6SHong Zhang }
851814e85d6SHong Zhang 
852814e85d6SHong Zhang /*@
853814e85d6SHong Zhang    TSGetCostGradients - Returns the gradients from the TSAdjointSolve()
854814e85d6SHong Zhang 
855814e85d6SHong Zhang    Not Collective, but Vec returned is parallel if TS is parallel
856814e85d6SHong Zhang 
857814e85d6SHong Zhang    Input Parameter:
858814e85d6SHong Zhang .  ts - the TS context obtained from TSCreate()
859814e85d6SHong Zhang 
860814e85d6SHong Zhang    Output Parameter:
861814e85d6SHong Zhang +  lambda - vectors containing the gradients of the cost functions with respect to the ODE/DAE solution variables
862814e85d6SHong Zhang -  mu - vectors containing the gradients of the cost functions with respect to the problem parameters
863814e85d6SHong Zhang 
864814e85d6SHong Zhang    Level: intermediate
865814e85d6SHong Zhang 
866b5b4017aSHong Zhang .seealso: TSSetCostGradients()
867814e85d6SHong Zhang @*/
868814e85d6SHong Zhang PetscErrorCode TSGetCostGradients(TS ts,PetscInt *numcost,Vec **lambda,Vec **mu)
869814e85d6SHong Zhang {
870814e85d6SHong Zhang   PetscFunctionBegin;
871814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
872814e85d6SHong Zhang   if (numcost) *numcost = ts->numcost;
873814e85d6SHong Zhang   if (lambda)  *lambda  = ts->vecs_sensi;
874814e85d6SHong Zhang   if (mu)      *mu      = ts->vecs_sensip;
875814e85d6SHong Zhang   PetscFunctionReturn(0);
876814e85d6SHong Zhang }
877814e85d6SHong Zhang 
878814e85d6SHong Zhang /*@
879b5b4017aSHong 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
880b5b4017aSHong Zhang       for use by the TSAdjoint routines.
881b5b4017aSHong Zhang 
882d083f849SBarry Smith    Logically Collective on TS
883b5b4017aSHong Zhang 
884b5b4017aSHong Zhang    Input Parameters:
885b5b4017aSHong Zhang +  ts - the TS context obtained from TSCreate()
886b5b4017aSHong Zhang .  numcost - number of cost functions
887b5b4017aSHong 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
888b5b4017aSHong 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
889b5b4017aSHong Zhang -  dir - the direction vector that are multiplied with the Hessian of the cost functions
890b5b4017aSHong Zhang 
891b5b4017aSHong Zhang    Level: beginner
892b5b4017aSHong Zhang 
893b5b4017aSHong Zhang    Notes: Hessian of the cost function is completely different from Hessian of the ODE/DAE system
894b5b4017aSHong Zhang 
895ecf68647SHong Zhang    For second-order adjoint, one needs to call this function and then TSAdjointSetForward() before TSSolve().
896b5b4017aSHong Zhang 
897b5b4017aSHong 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.
898b5b4017aSHong Zhang 
8993fca9621SHong Zhang    Passing NULL for lambda2 disables the second-order calculation.
900ecf68647SHong Zhang .seealso: TSAdjointSetForward()
901b5b4017aSHong Zhang @*/
902b5b4017aSHong Zhang PetscErrorCode TSSetCostHessianProducts(TS ts,PetscInt numcost,Vec *lambda2,Vec *mu2,Vec dir)
903b5b4017aSHong Zhang {
904b5b4017aSHong Zhang   PetscFunctionBegin;
905b5b4017aSHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
906b5b4017aSHong 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");
907b5b4017aSHong Zhang   ts->numcost       = numcost;
908b5b4017aSHong Zhang   ts->vecs_sensi2   = lambda2;
90933f72c5dSHong Zhang   ts->vecs_sensi2p  = mu2;
910b5b4017aSHong Zhang   ts->vec_dir       = dir;
911b5b4017aSHong Zhang   PetscFunctionReturn(0);
912b5b4017aSHong Zhang }
913b5b4017aSHong Zhang 
914b5b4017aSHong Zhang /*@
915b5b4017aSHong Zhang    TSGetCostHessianProducts - Returns the gradients from the TSAdjointSolve()
916b5b4017aSHong Zhang 
917b5b4017aSHong Zhang    Not Collective, but Vec returned is parallel if TS is parallel
918b5b4017aSHong Zhang 
919b5b4017aSHong Zhang    Input Parameter:
920b5b4017aSHong Zhang .  ts - the TS context obtained from TSCreate()
921b5b4017aSHong Zhang 
922b5b4017aSHong Zhang    Output Parameter:
923b5b4017aSHong Zhang +  numcost - number of cost functions
924b5b4017aSHong 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
925b5b4017aSHong 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
926b5b4017aSHong Zhang -  dir - the direction vector that are multiplied with the Hessian of the cost functions
927b5b4017aSHong Zhang 
928b5b4017aSHong Zhang    Level: intermediate
929b5b4017aSHong Zhang 
930b5b4017aSHong Zhang .seealso: TSSetCostHessianProducts()
931b5b4017aSHong Zhang @*/
932b5b4017aSHong Zhang PetscErrorCode TSGetCostHessianProducts(TS ts,PetscInt *numcost,Vec **lambda2,Vec **mu2, Vec *dir)
933b5b4017aSHong Zhang {
934b5b4017aSHong Zhang   PetscFunctionBegin;
935b5b4017aSHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
936b5b4017aSHong Zhang   if (numcost) *numcost = ts->numcost;
937b5b4017aSHong Zhang   if (lambda2) *lambda2 = ts->vecs_sensi2;
93833f72c5dSHong Zhang   if (mu2)     *mu2     = ts->vecs_sensi2p;
939b5b4017aSHong Zhang   if (dir)     *dir     = ts->vec_dir;
940b5b4017aSHong Zhang   PetscFunctionReturn(0);
941b5b4017aSHong Zhang }
942b5b4017aSHong Zhang 
943b5b4017aSHong Zhang /*@
944ecf68647SHong Zhang   TSAdjointSetForward - Trigger the tangent linear solver and initialize the forward sensitivities
945b5b4017aSHong Zhang 
946d083f849SBarry Smith   Logically Collective on TS
947b5b4017aSHong Zhang 
948b5b4017aSHong Zhang   Input Parameters:
949b5b4017aSHong Zhang +  ts - the TS context obtained from TSCreate()
950b5b4017aSHong Zhang -  didp - the derivative of initial values w.r.t. parameters
951b5b4017aSHong Zhang 
952b5b4017aSHong Zhang   Level: intermediate
953b5b4017aSHong Zhang 
954ecf68647SHong 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.
955b5b4017aSHong Zhang 
956ecf68647SHong Zhang .seealso: TSSetCostHessianProducts(), TSAdjointResetForward()
957b5b4017aSHong Zhang @*/
958ecf68647SHong Zhang PetscErrorCode TSAdjointSetForward(TS ts,Mat didp)
959b5b4017aSHong Zhang {
96033f72c5dSHong Zhang   Mat            A;
96133f72c5dSHong Zhang   Vec            sp;
96233f72c5dSHong Zhang   PetscScalar    *xarr;
96333f72c5dSHong Zhang   PetscInt       lsize;
964b5b4017aSHong Zhang   PetscErrorCode ierr;
965b5b4017aSHong Zhang 
966b5b4017aSHong Zhang   PetscFunctionBegin;
967b5b4017aSHong Zhang   ts->forward_solve = PETSC_TRUE; /* turn on tangent linear mode */
968b5b4017aSHong Zhang   if (!ts->vecs_sensi2) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must call TSSetCostHessianProducts() first");
96933f72c5dSHong Zhang   if (!ts->vec_dir) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Directional vector is missing. Call TSSetCostHessianProducts() to set it.");
97033f72c5dSHong Zhang   /* create a single-column dense matrix */
97133f72c5dSHong Zhang   ierr = VecGetLocalSize(ts->vec_sol,&lsize);CHKERRQ(ierr);
972f63bf25fSHong Zhang   ierr = MatCreateDense(PetscObjectComm((PetscObject)ts),lsize,PETSC_DECIDE,PETSC_DECIDE,1,NULL,&A);CHKERRQ(ierr);
97333f72c5dSHong Zhang 
97433f72c5dSHong Zhang   ierr = VecDuplicate(ts->vec_sol,&sp);CHKERRQ(ierr);
97533f72c5dSHong Zhang   ierr = MatDenseGetColumn(A,0,&xarr);CHKERRQ(ierr);
97633f72c5dSHong Zhang   ierr = VecPlaceArray(sp,xarr);CHKERRQ(ierr);
977ecf68647SHong Zhang   if (ts->vecs_sensi2p) { /* tangent linear variable initialized as 2*dIdP*dir */
97833f72c5dSHong Zhang     if (didp) {
97933f72c5dSHong Zhang       ierr = MatMult(didp,ts->vec_dir,sp);CHKERRQ(ierr);
98033f72c5dSHong Zhang       ierr = VecScale(sp,2.);CHKERRQ(ierr);
98133f72c5dSHong Zhang     } else {
98233f72c5dSHong Zhang       ierr = VecZeroEntries(sp);CHKERRQ(ierr);
98333f72c5dSHong Zhang     }
984ecf68647SHong Zhang   } else { /* tangent linear variable initialized as dir */
98533f72c5dSHong Zhang     ierr = VecCopy(ts->vec_dir,sp);CHKERRQ(ierr);
98633f72c5dSHong Zhang   }
98733f72c5dSHong Zhang   ierr = VecResetArray(sp);CHKERRQ(ierr);
98833f72c5dSHong Zhang   ierr = MatDenseRestoreColumn(A,&xarr);CHKERRQ(ierr);
98933f72c5dSHong Zhang   ierr = VecDestroy(&sp);CHKERRQ(ierr);
99033f72c5dSHong Zhang 
99133f72c5dSHong Zhang   ierr = TSForwardSetInitialSensitivities(ts,A);CHKERRQ(ierr); /* if didp is NULL, identity matrix is assumed */
99233f72c5dSHong Zhang 
99333f72c5dSHong Zhang   ierr = MatDestroy(&A);CHKERRQ(ierr);
994b5b4017aSHong Zhang   PetscFunctionReturn(0);
995b5b4017aSHong Zhang }
996b5b4017aSHong Zhang 
997b5b4017aSHong Zhang /*@
998ecf68647SHong Zhang   TSAdjointResetForward - Reset the tangent linear solver and destroy the tangent linear context
999ecf68647SHong Zhang 
1000d083f849SBarry Smith   Logically Collective on TS
1001ecf68647SHong Zhang 
1002ecf68647SHong Zhang   Input Parameters:
1003ecf68647SHong Zhang .  ts - the TS context obtained from TSCreate()
1004ecf68647SHong Zhang 
1005ecf68647SHong Zhang   Level: intermediate
1006ecf68647SHong Zhang 
1007ecf68647SHong Zhang .seealso: TSAdjointSetForward()
1008ecf68647SHong Zhang @*/
1009ecf68647SHong Zhang PetscErrorCode TSAdjointResetForward(TS ts)
1010ecf68647SHong Zhang {
1011ecf68647SHong Zhang   PetscErrorCode ierr;
1012ecf68647SHong Zhang 
1013ecf68647SHong Zhang   PetscFunctionBegin;
1014ecf68647SHong Zhang   ts->forward_solve = PETSC_FALSE; /* turn off tangent linear mode */
1015ecf68647SHong Zhang   ierr = TSForwardReset(ts);CHKERRQ(ierr);
1016ecf68647SHong Zhang   PetscFunctionReturn(0);
1017ecf68647SHong Zhang }
1018ecf68647SHong Zhang 
1019ecf68647SHong Zhang /*@
1020814e85d6SHong Zhang    TSAdjointSetUp - Sets up the internal data structures for the later use
1021814e85d6SHong Zhang    of an adjoint solver
1022814e85d6SHong Zhang 
1023814e85d6SHong Zhang    Collective on TS
1024814e85d6SHong Zhang 
1025814e85d6SHong Zhang    Input Parameter:
1026814e85d6SHong Zhang .  ts - the TS context obtained from TSCreate()
1027814e85d6SHong Zhang 
1028814e85d6SHong Zhang    Level: advanced
1029814e85d6SHong Zhang 
1030814e85d6SHong Zhang .seealso: TSCreate(), TSAdjointStep(), TSSetCostGradients()
1031814e85d6SHong Zhang @*/
1032814e85d6SHong Zhang PetscErrorCode TSAdjointSetUp(TS ts)
1033814e85d6SHong Zhang {
1034881c1a9bSHong Zhang   TSTrajectory     tj;
1035881c1a9bSHong Zhang   PetscBool        match;
1036814e85d6SHong Zhang   PetscErrorCode   ierr;
1037814e85d6SHong Zhang 
1038814e85d6SHong Zhang   PetscFunctionBegin;
1039814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1040814e85d6SHong Zhang   if (ts->adjointsetupcalled) PetscFunctionReturn(0);
1041814e85d6SHong Zhang   if (!ts->vecs_sensi) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetCostGradients() first");
104233f72c5dSHong Zhang   if (ts->vecs_sensip && !ts->Jacp && !ts->Jacprhs) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetRHSJacobianP() or TSSetIJacobianP() first");
1043881c1a9bSHong Zhang   ierr = TSGetTrajectory(ts,&tj);CHKERRQ(ierr);
10448e224c67SHong Zhang   ierr = PetscObjectTypeCompare((PetscObject)tj,TSTRAJECTORYBASIC,&match);CHKERRQ(ierr);
1045881c1a9bSHong Zhang   if (match) {
1046881c1a9bSHong Zhang     PetscBool solution_only;
1047881c1a9bSHong Zhang     ierr = TSTrajectoryGetSolutionOnly(tj,&solution_only);CHKERRQ(ierr);
1048881c1a9bSHong 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");
1049881c1a9bSHong Zhang   }
10509ffb3502SHong Zhang   ierr = TSTrajectorySetUseHistory(tj,PETSC_FALSE);CHKERRQ(ierr); /* not use TSHistory */
105133f72c5dSHong Zhang 
1052cd4cee2dSHong Zhang   if (ts->quadraturets) { /* if there is integral in the cost function */
1053cd4cee2dSHong Zhang     ierr = VecDuplicate(ts->vecs_sensi[0],&ts->vec_drdu_col);CHKERRQ(ierr);
1054814e85d6SHong Zhang     if (ts->vecs_sensip){
1055cd4cee2dSHong Zhang       ierr = VecDuplicate(ts->vecs_sensip[0],&ts->vec_drdp_col);CHKERRQ(ierr);
1056814e85d6SHong Zhang     }
1057814e85d6SHong Zhang   }
1058814e85d6SHong Zhang 
1059814e85d6SHong Zhang   if (ts->ops->adjointsetup) {
1060814e85d6SHong Zhang     ierr = (*ts->ops->adjointsetup)(ts);CHKERRQ(ierr);
1061814e85d6SHong Zhang   }
1062814e85d6SHong Zhang   ts->adjointsetupcalled = PETSC_TRUE;
1063814e85d6SHong Zhang   PetscFunctionReturn(0);
1064814e85d6SHong Zhang }
1065814e85d6SHong Zhang 
1066814e85d6SHong Zhang /*@
1067ece46509SHong Zhang    TSAdjointReset - Resets a TSAdjoint context and removes any allocated Vecs and Mats.
1068ece46509SHong Zhang 
1069ece46509SHong Zhang    Collective on TS
1070ece46509SHong Zhang 
1071ece46509SHong Zhang    Input Parameter:
1072ece46509SHong Zhang .  ts - the TS context obtained from TSCreate()
1073ece46509SHong Zhang 
1074ece46509SHong Zhang    Level: beginner
1075ece46509SHong Zhang 
10769ffb3502SHong Zhang .seealso: TSCreate(), TSAdjointSetUp(), TSADestroy()
1077ece46509SHong Zhang @*/
1078ece46509SHong Zhang PetscErrorCode TSAdjointReset(TS ts)
1079ece46509SHong Zhang {
1080ece46509SHong Zhang   PetscErrorCode ierr;
1081ece46509SHong Zhang 
1082ece46509SHong Zhang   PetscFunctionBegin;
1083ece46509SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1084ece46509SHong Zhang   if (ts->ops->adjointreset) {
1085ece46509SHong Zhang     ierr = (*ts->ops->adjointreset)(ts);CHKERRQ(ierr);
1086ece46509SHong Zhang   }
10877207e0fdSHong Zhang   if (ts->quadraturets) { /* if there is integral in the cost function */
10887207e0fdSHong Zhang     ierr = VecDestroy(&ts->vec_drdu_col);CHKERRQ(ierr);
10897207e0fdSHong Zhang     if (ts->vecs_sensip){
10907207e0fdSHong Zhang       ierr = VecDestroy(&ts->vec_drdp_col);CHKERRQ(ierr);
10917207e0fdSHong Zhang     }
10927207e0fdSHong Zhang   }
1093ece46509SHong Zhang   ts->vecs_sensi         = NULL;
1094ece46509SHong Zhang   ts->vecs_sensip        = NULL;
1095ece46509SHong Zhang   ts->vecs_sensi2        = NULL;
109633f72c5dSHong Zhang   ts->vecs_sensi2p       = NULL;
1097ece46509SHong Zhang   ts->vec_dir            = NULL;
1098ece46509SHong Zhang   ts->adjointsetupcalled = PETSC_FALSE;
1099ece46509SHong Zhang   PetscFunctionReturn(0);
1100ece46509SHong Zhang }
1101ece46509SHong Zhang 
1102ece46509SHong Zhang /*@
1103814e85d6SHong Zhang    TSAdjointSetSteps - Sets the number of steps the adjoint solver should take backward in time
1104814e85d6SHong Zhang 
1105814e85d6SHong Zhang    Logically Collective on TS
1106814e85d6SHong Zhang 
1107814e85d6SHong Zhang    Input Parameters:
1108814e85d6SHong Zhang +  ts - the TS context obtained from TSCreate()
1109a2b725a8SWilliam Gropp -  steps - number of steps to use
1110814e85d6SHong Zhang 
1111814e85d6SHong Zhang    Level: intermediate
1112814e85d6SHong Zhang 
1113814e85d6SHong Zhang    Notes:
1114814e85d6SHong Zhang     Normally one does not call this and TSAdjointSolve() integrates back to the original timestep. One can call this
1115814e85d6SHong Zhang           so as to integrate back to less than the original timestep
1116814e85d6SHong Zhang 
1117814e85d6SHong Zhang .seealso: TSSetExactFinalTime()
1118814e85d6SHong Zhang @*/
1119814e85d6SHong Zhang PetscErrorCode TSAdjointSetSteps(TS ts,PetscInt steps)
1120814e85d6SHong Zhang {
1121814e85d6SHong Zhang   PetscFunctionBegin;
1122814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1123814e85d6SHong Zhang   PetscValidLogicalCollectiveInt(ts,steps,2);
1124814e85d6SHong Zhang   if (steps < 0) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Cannot step back a negative number of steps");
1125814e85d6SHong Zhang   if (steps > ts->steps) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Cannot step back more than the total number of forward steps");
1126814e85d6SHong Zhang   ts->adjoint_max_steps = steps;
1127814e85d6SHong Zhang   PetscFunctionReturn(0);
1128814e85d6SHong Zhang }
1129814e85d6SHong Zhang 
1130814e85d6SHong Zhang /*@C
1131814e85d6SHong Zhang   TSAdjointSetRHSJacobian - Deprecated, use TSSetRHSJacobianP()
1132814e85d6SHong Zhang 
1133814e85d6SHong Zhang   Level: deprecated
1134814e85d6SHong Zhang 
1135814e85d6SHong Zhang @*/
1136814e85d6SHong Zhang PetscErrorCode TSAdjointSetRHSJacobian(TS ts,Mat Amat,PetscErrorCode (*func)(TS,PetscReal,Vec,Mat,void*),void *ctx)
1137814e85d6SHong Zhang {
1138814e85d6SHong Zhang   PetscErrorCode ierr;
1139814e85d6SHong Zhang 
1140814e85d6SHong Zhang   PetscFunctionBegin;
1141814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1142814e85d6SHong Zhang   PetscValidHeaderSpecific(Amat,MAT_CLASSID,2);
1143814e85d6SHong Zhang 
1144814e85d6SHong Zhang   ts->rhsjacobianp    = func;
1145814e85d6SHong Zhang   ts->rhsjacobianpctx = ctx;
1146814e85d6SHong Zhang   if(Amat) {
1147814e85d6SHong Zhang     ierr = PetscObjectReference((PetscObject)Amat);CHKERRQ(ierr);
1148814e85d6SHong Zhang     ierr = MatDestroy(&ts->Jacp);CHKERRQ(ierr);
1149814e85d6SHong Zhang     ts->Jacp = Amat;
1150814e85d6SHong Zhang   }
1151814e85d6SHong Zhang   PetscFunctionReturn(0);
1152814e85d6SHong Zhang }
1153814e85d6SHong Zhang 
1154814e85d6SHong Zhang /*@C
1155814e85d6SHong Zhang   TSAdjointComputeRHSJacobian - Deprecated, use TSComputeRHSJacobianP()
1156814e85d6SHong Zhang 
1157814e85d6SHong Zhang   Level: deprecated
1158814e85d6SHong Zhang 
1159814e85d6SHong Zhang @*/
1160c9ad9de0SHong Zhang PetscErrorCode TSAdjointComputeRHSJacobian(TS ts,PetscReal t,Vec U,Mat Amat)
1161814e85d6SHong Zhang {
1162814e85d6SHong Zhang   PetscErrorCode ierr;
1163814e85d6SHong Zhang 
1164814e85d6SHong Zhang   PetscFunctionBegin;
1165814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1166c9ad9de0SHong Zhang   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
1167814e85d6SHong Zhang   PetscValidPointer(Amat,4);
1168814e85d6SHong Zhang 
1169814e85d6SHong Zhang   PetscStackPush("TS user JacobianP function for sensitivity analysis");
1170c9ad9de0SHong Zhang   ierr = (*ts->rhsjacobianp)(ts,t,U,Amat,ts->rhsjacobianpctx); CHKERRQ(ierr);
1171814e85d6SHong Zhang   PetscStackPop;
1172814e85d6SHong Zhang   PetscFunctionReturn(0);
1173814e85d6SHong Zhang }
1174814e85d6SHong Zhang 
1175814e85d6SHong Zhang /*@
1176d76be551SHong Zhang   TSAdjointComputeDRDYFunction - Deprecated, use TSGetQuadratureTS() then TSComputeRHSJacobian()
1177814e85d6SHong Zhang 
1178814e85d6SHong Zhang   Level: deprecated
1179814e85d6SHong Zhang 
1180814e85d6SHong Zhang @*/
1181c9ad9de0SHong Zhang PetscErrorCode TSAdjointComputeDRDYFunction(TS ts,PetscReal t,Vec U,Vec *DRDU)
1182814e85d6SHong Zhang {
1183814e85d6SHong Zhang   PetscErrorCode ierr;
1184814e85d6SHong Zhang 
1185814e85d6SHong Zhang   PetscFunctionBegin;
1186814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1187c9ad9de0SHong Zhang   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
1188814e85d6SHong Zhang 
1189814e85d6SHong Zhang   PetscStackPush("TS user DRDY function for sensitivity analysis");
1190c9ad9de0SHong Zhang   ierr = (*ts->drdufunction)(ts,t,U,DRDU,ts->costintegrandctx); CHKERRQ(ierr);
1191814e85d6SHong Zhang   PetscStackPop;
1192814e85d6SHong Zhang   PetscFunctionReturn(0);
1193814e85d6SHong Zhang }
1194814e85d6SHong Zhang 
1195814e85d6SHong Zhang /*@
1196d76be551SHong Zhang   TSAdjointComputeDRDPFunction - Deprecated, use TSGetQuadratureTS() then TSComputeRHSJacobianP()
1197814e85d6SHong Zhang 
1198814e85d6SHong Zhang   Level: deprecated
1199814e85d6SHong Zhang 
1200814e85d6SHong Zhang @*/
1201c9ad9de0SHong Zhang PetscErrorCode TSAdjointComputeDRDPFunction(TS ts,PetscReal t,Vec U,Vec *DRDP)
1202814e85d6SHong Zhang {
1203814e85d6SHong Zhang   PetscErrorCode ierr;
1204814e85d6SHong Zhang 
1205814e85d6SHong Zhang   PetscFunctionBegin;
1206814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1207c9ad9de0SHong Zhang   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
1208814e85d6SHong Zhang 
1209814e85d6SHong Zhang   PetscStackPush("TS user DRDP function for sensitivity analysis");
1210c9ad9de0SHong Zhang   ierr = (*ts->drdpfunction)(ts,t,U,DRDP,ts->costintegrandctx); CHKERRQ(ierr);
1211814e85d6SHong Zhang   PetscStackPop;
1212814e85d6SHong Zhang   PetscFunctionReturn(0);
1213814e85d6SHong Zhang }
1214814e85d6SHong Zhang 
1215814e85d6SHong Zhang /*@C
1216814e85d6SHong Zhang    TSAdjointMonitorSensi - monitors the first lambda sensitivity
1217814e85d6SHong Zhang 
1218814e85d6SHong Zhang    Level: intermediate
1219814e85d6SHong Zhang 
1220814e85d6SHong Zhang .seealso: TSAdjointMonitorSet()
1221814e85d6SHong Zhang @*/
1222814e85d6SHong Zhang PetscErrorCode TSAdjointMonitorSensi(TS ts,PetscInt step,PetscReal ptime,Vec v,PetscInt numcost,Vec *lambda,Vec *mu,PetscViewerAndFormat *vf)
1223814e85d6SHong Zhang {
1224814e85d6SHong Zhang   PetscErrorCode ierr;
1225814e85d6SHong Zhang   PetscViewer    viewer = vf->viewer;
1226814e85d6SHong Zhang 
1227814e85d6SHong Zhang   PetscFunctionBegin;
1228814e85d6SHong Zhang   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4);
1229814e85d6SHong Zhang   ierr = PetscViewerPushFormat(viewer,vf->format);CHKERRQ(ierr);
1230814e85d6SHong Zhang   ierr = VecView(lambda[0],viewer);CHKERRQ(ierr);
1231814e85d6SHong Zhang   ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
1232814e85d6SHong Zhang   PetscFunctionReturn(0);
1233814e85d6SHong Zhang }
1234814e85d6SHong Zhang 
1235814e85d6SHong Zhang /*@C
1236814e85d6SHong Zhang    TSAdjointMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user
1237814e85d6SHong Zhang 
1238814e85d6SHong Zhang    Collective on TS
1239814e85d6SHong Zhang 
1240814e85d6SHong Zhang    Input Parameters:
1241814e85d6SHong Zhang +  ts - TS object you wish to monitor
1242814e85d6SHong Zhang .  name - the monitor type one is seeking
1243814e85d6SHong Zhang .  help - message indicating what monitoring is done
1244814e85d6SHong Zhang .  manual - manual page for the monitor
1245814e85d6SHong Zhang .  monitor - the monitor function
1246814e85d6SHong 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
1247814e85d6SHong Zhang 
1248814e85d6SHong Zhang    Level: developer
1249814e85d6SHong Zhang 
1250814e85d6SHong Zhang .seealso: PetscOptionsGetViewer(), PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(),
1251814e85d6SHong Zhang           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
1252814e85d6SHong Zhang           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
1253814e85d6SHong Zhang           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
1254814e85d6SHong Zhang           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1255814e85d6SHong Zhang           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1256814e85d6SHong Zhang           PetscOptionsFList(), PetscOptionsEList()
1257814e85d6SHong Zhang @*/
1258814e85d6SHong 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*))
1259814e85d6SHong Zhang {
1260814e85d6SHong Zhang   PetscErrorCode    ierr;
1261814e85d6SHong Zhang   PetscViewer       viewer;
1262814e85d6SHong Zhang   PetscViewerFormat format;
1263814e85d6SHong Zhang   PetscBool         flg;
1264814e85d6SHong Zhang 
1265814e85d6SHong Zhang   PetscFunctionBegin;
126616413a6aSBarry Smith   ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)ts),((PetscObject) ts)->options,((PetscObject)ts)->prefix,name,&viewer,&format,&flg);CHKERRQ(ierr);
1267814e85d6SHong Zhang   if (flg) {
1268814e85d6SHong Zhang     PetscViewerAndFormat *vf;
1269814e85d6SHong Zhang     ierr = PetscViewerAndFormatCreate(viewer,format,&vf);CHKERRQ(ierr);
1270814e85d6SHong Zhang     ierr = PetscObjectDereference((PetscObject)viewer);CHKERRQ(ierr);
1271814e85d6SHong Zhang     if (monitorsetup) {
1272814e85d6SHong Zhang       ierr = (*monitorsetup)(ts,vf);CHKERRQ(ierr);
1273814e85d6SHong Zhang     }
1274814e85d6SHong Zhang     ierr = TSAdjointMonitorSet(ts,(PetscErrorCode (*)(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*))monitor,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy);CHKERRQ(ierr);
1275814e85d6SHong Zhang   }
1276814e85d6SHong Zhang   PetscFunctionReturn(0);
1277814e85d6SHong Zhang }
1278814e85d6SHong Zhang 
1279814e85d6SHong Zhang /*@C
1280814e85d6SHong Zhang    TSAdjointMonitorSet - Sets an ADDITIONAL function that is to be used at every
1281814e85d6SHong Zhang    timestep to display the iteration's  progress.
1282814e85d6SHong Zhang 
1283814e85d6SHong Zhang    Logically Collective on TS
1284814e85d6SHong Zhang 
1285814e85d6SHong Zhang    Input Parameters:
1286814e85d6SHong Zhang +  ts - the TS context obtained from TSCreate()
1287814e85d6SHong Zhang .  adjointmonitor - monitoring routine
1288814e85d6SHong Zhang .  adjointmctx - [optional] user-defined context for private data for the
1289814e85d6SHong Zhang              monitor routine (use NULL if no context is desired)
1290814e85d6SHong Zhang -  adjointmonitordestroy - [optional] routine that frees monitor context
1291814e85d6SHong Zhang           (may be NULL)
1292814e85d6SHong Zhang 
1293814e85d6SHong Zhang    Calling sequence of monitor:
1294814e85d6SHong Zhang $    int adjointmonitor(TS ts,PetscInt steps,PetscReal time,Vec u,PetscInt numcost,Vec *lambda, Vec *mu,void *adjointmctx)
1295814e85d6SHong Zhang 
1296814e85d6SHong Zhang +    ts - the TS context
1297814e85d6SHong 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
1298814e85d6SHong Zhang                                been interpolated to)
1299814e85d6SHong Zhang .    time - current time
1300814e85d6SHong Zhang .    u - current iterate
1301814e85d6SHong Zhang .    numcost - number of cost functionos
1302814e85d6SHong Zhang .    lambda - sensitivities to initial conditions
1303814e85d6SHong Zhang .    mu - sensitivities to parameters
1304814e85d6SHong Zhang -    adjointmctx - [optional] adjoint monitoring context
1305814e85d6SHong Zhang 
1306814e85d6SHong Zhang    Notes:
1307814e85d6SHong Zhang    This routine adds an additional monitor to the list of monitors that
1308814e85d6SHong Zhang    already has been loaded.
1309814e85d6SHong Zhang 
1310814e85d6SHong Zhang    Fortran Notes:
1311814e85d6SHong Zhang     Only a single monitor function can be set for each TS object
1312814e85d6SHong Zhang 
1313814e85d6SHong Zhang    Level: intermediate
1314814e85d6SHong Zhang 
1315814e85d6SHong Zhang .seealso: TSAdjointMonitorCancel()
1316814e85d6SHong Zhang @*/
1317814e85d6SHong Zhang PetscErrorCode TSAdjointMonitorSet(TS ts,PetscErrorCode (*adjointmonitor)(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*),void *adjointmctx,PetscErrorCode (*adjointmdestroy)(void**))
1318814e85d6SHong Zhang {
1319814e85d6SHong Zhang   PetscErrorCode ierr;
1320814e85d6SHong Zhang   PetscInt       i;
1321814e85d6SHong Zhang   PetscBool      identical;
1322814e85d6SHong Zhang 
1323814e85d6SHong Zhang   PetscFunctionBegin;
1324814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1325814e85d6SHong Zhang   for (i=0; i<ts->numbermonitors;i++) {
1326814e85d6SHong Zhang     ierr = PetscMonitorCompare((PetscErrorCode (*)(void))adjointmonitor,adjointmctx,adjointmdestroy,(PetscErrorCode (*)(void))ts->adjointmonitor[i],ts->adjointmonitorcontext[i],ts->adjointmonitordestroy[i],&identical);CHKERRQ(ierr);
1327814e85d6SHong Zhang     if (identical) PetscFunctionReturn(0);
1328814e85d6SHong Zhang   }
1329814e85d6SHong Zhang   if (ts->numberadjointmonitors >= MAXTSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many adjoint monitors set");
1330814e85d6SHong Zhang   ts->adjointmonitor[ts->numberadjointmonitors]          = adjointmonitor;
1331814e85d6SHong Zhang   ts->adjointmonitordestroy[ts->numberadjointmonitors]   = adjointmdestroy;
1332814e85d6SHong Zhang   ts->adjointmonitorcontext[ts->numberadjointmonitors++] = (void*)adjointmctx;
1333814e85d6SHong Zhang   PetscFunctionReturn(0);
1334814e85d6SHong Zhang }
1335814e85d6SHong Zhang 
1336814e85d6SHong Zhang /*@C
1337814e85d6SHong Zhang    TSAdjointMonitorCancel - Clears all the adjoint monitors that have been set on a time-step object.
1338814e85d6SHong Zhang 
1339814e85d6SHong Zhang    Logically Collective on TS
1340814e85d6SHong Zhang 
1341814e85d6SHong Zhang    Input Parameters:
1342814e85d6SHong Zhang .  ts - the TS context obtained from TSCreate()
1343814e85d6SHong Zhang 
1344814e85d6SHong Zhang    Notes:
1345814e85d6SHong Zhang    There is no way to remove a single, specific monitor.
1346814e85d6SHong Zhang 
1347814e85d6SHong Zhang    Level: intermediate
1348814e85d6SHong Zhang 
1349814e85d6SHong Zhang .seealso: TSAdjointMonitorSet()
1350814e85d6SHong Zhang @*/
1351814e85d6SHong Zhang PetscErrorCode TSAdjointMonitorCancel(TS ts)
1352814e85d6SHong Zhang {
1353814e85d6SHong Zhang   PetscErrorCode ierr;
1354814e85d6SHong Zhang   PetscInt       i;
1355814e85d6SHong Zhang 
1356814e85d6SHong Zhang   PetscFunctionBegin;
1357814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1358814e85d6SHong Zhang   for (i=0; i<ts->numberadjointmonitors; i++) {
1359814e85d6SHong Zhang     if (ts->adjointmonitordestroy[i]) {
1360814e85d6SHong Zhang       ierr = (*ts->adjointmonitordestroy[i])(&ts->adjointmonitorcontext[i]);CHKERRQ(ierr);
1361814e85d6SHong Zhang     }
1362814e85d6SHong Zhang   }
1363814e85d6SHong Zhang   ts->numberadjointmonitors = 0;
1364814e85d6SHong Zhang   PetscFunctionReturn(0);
1365814e85d6SHong Zhang }
1366814e85d6SHong Zhang 
1367814e85d6SHong Zhang /*@C
1368814e85d6SHong Zhang    TSAdjointMonitorDefault - the default monitor of adjoint computations
1369814e85d6SHong Zhang 
1370814e85d6SHong Zhang    Level: intermediate
1371814e85d6SHong Zhang 
1372814e85d6SHong Zhang .seealso: TSAdjointMonitorSet()
1373814e85d6SHong Zhang @*/
1374814e85d6SHong Zhang PetscErrorCode TSAdjointMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,PetscInt numcost,Vec *lambda,Vec *mu,PetscViewerAndFormat *vf)
1375814e85d6SHong Zhang {
1376814e85d6SHong Zhang   PetscErrorCode ierr;
1377814e85d6SHong Zhang   PetscViewer    viewer = vf->viewer;
1378814e85d6SHong Zhang 
1379814e85d6SHong Zhang   PetscFunctionBegin;
1380814e85d6SHong Zhang   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4);
1381814e85d6SHong Zhang   ierr = PetscViewerPushFormat(viewer,vf->format);CHKERRQ(ierr);
1382814e85d6SHong Zhang   ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
1383814e85d6SHong 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);
1384814e85d6SHong Zhang   ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
1385814e85d6SHong Zhang   ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
1386814e85d6SHong Zhang   PetscFunctionReturn(0);
1387814e85d6SHong Zhang }
1388814e85d6SHong Zhang 
1389814e85d6SHong Zhang /*@C
1390814e85d6SHong Zhang    TSAdjointMonitorDrawSensi - Monitors progress of the adjoint TS solvers by calling
1391814e85d6SHong Zhang    VecView() for the sensitivities to initial states at each timestep
1392814e85d6SHong Zhang 
1393814e85d6SHong Zhang    Collective on TS
1394814e85d6SHong Zhang 
1395814e85d6SHong Zhang    Input Parameters:
1396814e85d6SHong Zhang +  ts - the TS context
1397814e85d6SHong Zhang .  step - current time-step
1398814e85d6SHong Zhang .  ptime - current time
1399814e85d6SHong Zhang .  u - current state
1400814e85d6SHong Zhang .  numcost - number of cost functions
1401814e85d6SHong Zhang .  lambda - sensitivities to initial conditions
1402814e85d6SHong Zhang .  mu - sensitivities to parameters
1403814e85d6SHong Zhang -  dummy - either a viewer or NULL
1404814e85d6SHong Zhang 
1405814e85d6SHong Zhang    Level: intermediate
1406814e85d6SHong Zhang 
1407814e85d6SHong Zhang .seealso: TSAdjointMonitorSet(), TSAdjointMonitorDefault(), VecView()
1408814e85d6SHong Zhang @*/
1409814e85d6SHong Zhang PetscErrorCode TSAdjointMonitorDrawSensi(TS ts,PetscInt step,PetscReal ptime,Vec u,PetscInt numcost,Vec *lambda,Vec *mu,void *dummy)
1410814e85d6SHong Zhang {
1411814e85d6SHong Zhang   PetscErrorCode   ierr;
1412814e85d6SHong Zhang   TSMonitorDrawCtx ictx = (TSMonitorDrawCtx)dummy;
1413814e85d6SHong Zhang   PetscDraw        draw;
1414814e85d6SHong Zhang   PetscReal        xl,yl,xr,yr,h;
1415814e85d6SHong Zhang   char             time[32];
1416814e85d6SHong Zhang 
1417814e85d6SHong Zhang   PetscFunctionBegin;
1418814e85d6SHong Zhang   if (!(((ictx->howoften > 0) && (!(step % ictx->howoften))) || ((ictx->howoften == -1) && ts->reason))) PetscFunctionReturn(0);
1419814e85d6SHong Zhang 
1420814e85d6SHong Zhang   ierr = VecView(lambda[0],ictx->viewer);CHKERRQ(ierr);
1421814e85d6SHong Zhang   ierr = PetscViewerDrawGetDraw(ictx->viewer,0,&draw);CHKERRQ(ierr);
1422814e85d6SHong Zhang   ierr = PetscSNPrintf(time,32,"Timestep %d Time %g",(int)step,(double)ptime);CHKERRQ(ierr);
1423814e85d6SHong Zhang   ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
1424814e85d6SHong Zhang   h    = yl + .95*(yr - yl);
1425814e85d6SHong Zhang   ierr = PetscDrawStringCentered(draw,.5*(xl+xr),h,PETSC_DRAW_BLACK,time);CHKERRQ(ierr);
1426814e85d6SHong Zhang   ierr = PetscDrawFlush(draw);CHKERRQ(ierr);
1427814e85d6SHong Zhang   PetscFunctionReturn(0);
1428814e85d6SHong Zhang }
1429814e85d6SHong Zhang 
1430814e85d6SHong Zhang /*
1431814e85d6SHong Zhang    TSAdjointSetFromOptions - Sets various TSAdjoint parameters from user options.
1432814e85d6SHong Zhang 
1433814e85d6SHong Zhang    Collective on TSAdjoint
1434814e85d6SHong Zhang 
1435814e85d6SHong Zhang    Input Parameter:
1436814e85d6SHong Zhang .  ts - the TS context
1437814e85d6SHong Zhang 
1438814e85d6SHong Zhang    Options Database Keys:
1439814e85d6SHong Zhang +  -ts_adjoint_solve <yes,no> After solving the ODE/DAE solve the adjoint problem (requires -ts_save_trajectory)
1440814e85d6SHong Zhang .  -ts_adjoint_monitor - print information at each adjoint time step
1441814e85d6SHong Zhang -  -ts_adjoint_monitor_draw_sensi - monitor the sensitivity of the first cost function wrt initial conditions (lambda[0]) graphically
1442814e85d6SHong Zhang 
1443814e85d6SHong Zhang    Level: developer
1444814e85d6SHong Zhang 
1445814e85d6SHong Zhang    Notes:
1446814e85d6SHong Zhang     This is not normally called directly by users
1447814e85d6SHong Zhang 
1448814e85d6SHong Zhang .seealso: TSSetSaveTrajectory(), TSTrajectorySetUp()
1449814e85d6SHong Zhang */
1450814e85d6SHong Zhang PetscErrorCode TSAdjointSetFromOptions(PetscOptionItems *PetscOptionsObject,TS ts)
1451814e85d6SHong Zhang {
1452814e85d6SHong Zhang   PetscBool      tflg,opt;
1453814e85d6SHong Zhang   PetscErrorCode ierr;
1454814e85d6SHong Zhang 
1455814e85d6SHong Zhang   PetscFunctionBegin;
1456814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,2);
1457814e85d6SHong Zhang   ierr = PetscOptionsHead(PetscOptionsObject,"TS Adjoint options");CHKERRQ(ierr);
1458814e85d6SHong Zhang   tflg = ts->adjoint_solve ? PETSC_TRUE : PETSC_FALSE;
1459814e85d6SHong Zhang   ierr = PetscOptionsBool("-ts_adjoint_solve","Solve the adjoint problem immediately after solving the forward problem","",tflg,&tflg,&opt);CHKERRQ(ierr);
1460814e85d6SHong Zhang   if (opt) {
1461814e85d6SHong Zhang     ierr = TSSetSaveTrajectory(ts);CHKERRQ(ierr);
1462814e85d6SHong Zhang     ts->adjoint_solve = tflg;
1463814e85d6SHong Zhang   }
1464814e85d6SHong Zhang   ierr = TSAdjointMonitorSetFromOptions(ts,"-ts_adjoint_monitor","Monitor adjoint timestep size","TSAdjointMonitorDefault",TSAdjointMonitorDefault,NULL);CHKERRQ(ierr);
1465814e85d6SHong Zhang   ierr = TSAdjointMonitorSetFromOptions(ts,"-ts_adjoint_monitor_sensi","Monitor sensitivity in the adjoint computation","TSAdjointMonitorSensi",TSAdjointMonitorSensi,NULL);CHKERRQ(ierr);
1466814e85d6SHong Zhang   opt  = PETSC_FALSE;
1467814e85d6SHong Zhang   ierr = PetscOptionsName("-ts_adjoint_monitor_draw_sensi","Monitor adjoint sensitivities (lambda only) graphically","TSAdjointMonitorDrawSensi",&opt);CHKERRQ(ierr);
1468814e85d6SHong Zhang   if (opt) {
1469814e85d6SHong Zhang     TSMonitorDrawCtx ctx;
1470814e85d6SHong Zhang     PetscInt         howoften = 1;
1471814e85d6SHong Zhang 
1472814e85d6SHong Zhang     ierr = PetscOptionsInt("-ts_adjoint_monitor_draw_sensi","Monitor adjoint sensitivities (lambda only) graphically","TSAdjointMonitorDrawSensi",howoften,&howoften,NULL);CHKERRQ(ierr);
1473814e85d6SHong Zhang     ierr = TSMonitorDrawCtxCreate(PetscObjectComm((PetscObject)ts),0,0,PETSC_DECIDE,PETSC_DECIDE,300,300,howoften,&ctx);CHKERRQ(ierr);
1474814e85d6SHong Zhang     ierr = TSAdjointMonitorSet(ts,TSAdjointMonitorDrawSensi,ctx,(PetscErrorCode (*)(void**))TSMonitorDrawCtxDestroy);CHKERRQ(ierr);
1475814e85d6SHong Zhang   }
1476814e85d6SHong Zhang   PetscFunctionReturn(0);
1477814e85d6SHong Zhang }
1478814e85d6SHong Zhang 
1479814e85d6SHong Zhang /*@
1480814e85d6SHong Zhang    TSAdjointStep - Steps one time step backward in the adjoint run
1481814e85d6SHong Zhang 
1482814e85d6SHong Zhang    Collective on TS
1483814e85d6SHong Zhang 
1484814e85d6SHong Zhang    Input Parameter:
1485814e85d6SHong Zhang .  ts - the TS context obtained from TSCreate()
1486814e85d6SHong Zhang 
1487814e85d6SHong Zhang    Level: intermediate
1488814e85d6SHong Zhang 
1489814e85d6SHong Zhang .seealso: TSAdjointSetUp(), TSAdjointSolve()
1490814e85d6SHong Zhang @*/
1491814e85d6SHong Zhang PetscErrorCode TSAdjointStep(TS ts)
1492814e85d6SHong Zhang {
1493814e85d6SHong Zhang   DM               dm;
1494814e85d6SHong Zhang   PetscErrorCode   ierr;
1495814e85d6SHong Zhang 
1496814e85d6SHong Zhang   PetscFunctionBegin;
1497814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1498814e85d6SHong Zhang   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
1499814e85d6SHong Zhang   ierr = TSAdjointSetUp(ts);CHKERRQ(ierr);
1500814e85d6SHong Zhang 
1501814e85d6SHong Zhang   ts->reason = TS_CONVERGED_ITERATING;
1502814e85d6SHong Zhang   ts->ptime_prev = ts->ptime;
1503814e85d6SHong 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);
1504814e85d6SHong Zhang   ierr = PetscLogEventBegin(TS_AdjointStep,ts,0,0,0);CHKERRQ(ierr);
1505814e85d6SHong Zhang   ierr = (*ts->ops->adjointstep)(ts);CHKERRQ(ierr);
1506814e85d6SHong Zhang   ierr = PetscLogEventEnd(TS_AdjointStep,ts,0,0,0);CHKERRQ(ierr);
1507814e85d6SHong Zhang   ts->adjoint_steps++; ts->steps--;
1508814e85d6SHong Zhang 
1509814e85d6SHong Zhang   if (ts->reason < 0) {
1510814e85d6SHong Zhang     if (ts->errorifstepfailed) {
151105c9950eSHong 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]);
151205c9950eSHong Zhang       else SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_NOT_CONVERGED,"TSAdjointStep has failed due to %s",TSConvergedReasons[ts->reason]);
1513814e85d6SHong Zhang     }
1514814e85d6SHong Zhang   } else if (!ts->reason) {
1515814e85d6SHong Zhang     if (ts->adjoint_steps >= ts->adjoint_max_steps) ts->reason = TS_CONVERGED_ITS;
1516814e85d6SHong Zhang   }
1517814e85d6SHong Zhang   PetscFunctionReturn(0);
1518814e85d6SHong Zhang }
1519814e85d6SHong Zhang 
1520814e85d6SHong Zhang /*@
1521814e85d6SHong Zhang    TSAdjointSolve - Solves the discrete ajoint problem for an ODE/DAE
1522814e85d6SHong Zhang 
1523814e85d6SHong Zhang    Collective on TS
1524814e85d6SHong Zhang 
1525814e85d6SHong Zhang    Input Parameter:
1526814e85d6SHong Zhang .  ts - the TS context obtained from TSCreate()
1527814e85d6SHong Zhang 
1528814e85d6SHong Zhang    Options Database:
1529814e85d6SHong Zhang . -ts_adjoint_view_solution <viewerinfo> - views the first gradient with respect to the initial values
1530814e85d6SHong Zhang 
1531814e85d6SHong Zhang    Level: intermediate
1532814e85d6SHong Zhang 
1533814e85d6SHong Zhang    Notes:
1534814e85d6SHong Zhang    This must be called after a call to TSSolve() that solves the forward problem
1535814e85d6SHong Zhang 
1536814e85d6SHong Zhang    By default this will integrate back to the initial time, one can use TSAdjointSetSteps() to step back to a later time
1537814e85d6SHong Zhang 
1538814e85d6SHong Zhang .seealso: TSCreate(), TSSetCostGradients(), TSSetSolution(), TSAdjointStep()
1539814e85d6SHong Zhang @*/
1540814e85d6SHong Zhang PetscErrorCode TSAdjointSolve(TS ts)
1541814e85d6SHong Zhang {
15427f59fb53SHong Zhang #if defined(TSADJOINT_STAGE)
15437f59fb53SHong Zhang   PetscLogStage  adjoint_stage;
15447f59fb53SHong Zhang #endif
1545814e85d6SHong Zhang   PetscErrorCode ierr;
1546814e85d6SHong Zhang 
1547814e85d6SHong Zhang   PetscFunctionBegin;
1548814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
15497f59fb53SHong Zhang #if defined(TSADJOINT_STAGE)
15507f59fb53SHong Zhang   ierr = PetscLogStageRegister("TSAdjoint",&adjoint_stage);CHKERRQ(ierr);
15517f59fb53SHong Zhang   ierr = PetscLogStagePush(adjoint_stage);CHKERRQ(ierr);
15527f59fb53SHong Zhang #endif
1553814e85d6SHong Zhang   ierr = TSAdjointSetUp(ts);CHKERRQ(ierr);
1554814e85d6SHong Zhang 
1555814e85d6SHong Zhang   /* reset time step and iteration counters */
1556814e85d6SHong Zhang   ts->adjoint_steps     = 0;
1557814e85d6SHong Zhang   ts->ksp_its           = 0;
1558814e85d6SHong Zhang   ts->snes_its          = 0;
1559814e85d6SHong Zhang   ts->num_snes_failures = 0;
1560814e85d6SHong Zhang   ts->reject            = 0;
1561814e85d6SHong Zhang   ts->reason            = TS_CONVERGED_ITERATING;
1562814e85d6SHong Zhang 
1563814e85d6SHong Zhang   if (!ts->adjoint_max_steps) ts->adjoint_max_steps = ts->steps;
1564814e85d6SHong Zhang   if (ts->adjoint_steps >= ts->adjoint_max_steps) ts->reason = TS_CONVERGED_ITS;
1565814e85d6SHong Zhang 
1566814e85d6SHong Zhang   while (!ts->reason) {
1567814e85d6SHong Zhang     ierr = TSTrajectoryGet(ts->trajectory,ts,ts->steps,&ts->ptime);CHKERRQ(ierr);
1568814e85d6SHong Zhang     ierr = TSAdjointMonitor(ts,ts->steps,ts->ptime,ts->vec_sol,ts->numcost,ts->vecs_sensi,ts->vecs_sensip);CHKERRQ(ierr);
1569814e85d6SHong Zhang     ierr = TSAdjointEventHandler(ts);CHKERRQ(ierr);
1570814e85d6SHong Zhang     ierr = TSAdjointStep(ts);CHKERRQ(ierr);
1571cd4cee2dSHong Zhang     if ((ts->vec_costintegral || ts->quadraturets) && !ts->costintegralfwd) {
1572814e85d6SHong Zhang       ierr = TSAdjointCostIntegral(ts);CHKERRQ(ierr);
1573814e85d6SHong Zhang     }
1574814e85d6SHong Zhang   }
1575814e85d6SHong Zhang   ierr = TSTrajectoryGet(ts->trajectory,ts,ts->steps,&ts->ptime);CHKERRQ(ierr);
1576814e85d6SHong Zhang   ierr = TSAdjointMonitor(ts,ts->steps,ts->ptime,ts->vec_sol,ts->numcost,ts->vecs_sensi,ts->vecs_sensip);CHKERRQ(ierr);
1577814e85d6SHong Zhang   ts->solvetime = ts->ptime;
1578814e85d6SHong Zhang   ierr = TSTrajectoryViewFromOptions(ts->trajectory,NULL,"-ts_trajectory_view");CHKERRQ(ierr);
1579814e85d6SHong Zhang   ierr = VecViewFromOptions(ts->vecs_sensi[0],(PetscObject) ts, "-ts_adjoint_view_solution");CHKERRQ(ierr);
1580814e85d6SHong Zhang   ts->adjoint_max_steps = 0;
15817f59fb53SHong Zhang #if defined(TSADJOINT_STAGE)
15827f59fb53SHong Zhang   ierr = PetscLogStagePop();CHKERRQ(ierr);
15837f59fb53SHong Zhang #endif
1584814e85d6SHong Zhang   PetscFunctionReturn(0);
1585814e85d6SHong Zhang }
1586814e85d6SHong Zhang 
1587814e85d6SHong Zhang /*@C
1588814e85d6SHong Zhang    TSAdjointMonitor - Runs all user-provided adjoint monitor routines set using TSAdjointMonitorSet()
1589814e85d6SHong Zhang 
1590814e85d6SHong Zhang    Collective on TS
1591814e85d6SHong Zhang 
1592814e85d6SHong Zhang    Input Parameters:
1593814e85d6SHong Zhang +  ts - time stepping context obtained from TSCreate()
1594814e85d6SHong Zhang .  step - step number that has just completed
1595814e85d6SHong Zhang .  ptime - model time of the state
1596814e85d6SHong Zhang .  u - state at the current model time
1597814e85d6SHong Zhang .  numcost - number of cost functions (dimension of lambda  or mu)
1598814e85d6SHong Zhang .  lambda - vectors containing the gradients of the cost functions with respect to the ODE/DAE solution variables
1599814e85d6SHong Zhang -  mu - vectors containing the gradients of the cost functions with respect to the problem parameters
1600814e85d6SHong Zhang 
1601814e85d6SHong Zhang    Notes:
1602814e85d6SHong Zhang    TSAdjointMonitor() is typically used automatically within the time stepping implementations.
1603814e85d6SHong Zhang    Users would almost never call this routine directly.
1604814e85d6SHong Zhang 
1605814e85d6SHong Zhang    Level: developer
1606814e85d6SHong Zhang 
1607814e85d6SHong Zhang @*/
1608814e85d6SHong Zhang PetscErrorCode TSAdjointMonitor(TS ts,PetscInt step,PetscReal ptime,Vec u,PetscInt numcost,Vec *lambda, Vec *mu)
1609814e85d6SHong Zhang {
1610814e85d6SHong Zhang   PetscErrorCode ierr;
1611814e85d6SHong Zhang   PetscInt       i,n = ts->numberadjointmonitors;
1612814e85d6SHong Zhang 
1613814e85d6SHong Zhang   PetscFunctionBegin;
1614814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1615814e85d6SHong Zhang   PetscValidHeaderSpecific(u,VEC_CLASSID,4);
16168860a134SJunchao Zhang   ierr = VecLockReadPush(u);CHKERRQ(ierr);
1617814e85d6SHong Zhang   for (i=0; i<n; i++) {
1618814e85d6SHong Zhang     ierr = (*ts->adjointmonitor[i])(ts,step,ptime,u,numcost,lambda,mu,ts->adjointmonitorcontext[i]);CHKERRQ(ierr);
1619814e85d6SHong Zhang   }
16208860a134SJunchao Zhang   ierr = VecLockReadPop(u);CHKERRQ(ierr);
1621814e85d6SHong Zhang   PetscFunctionReturn(0);
1622814e85d6SHong Zhang }
1623814e85d6SHong Zhang 
1624814e85d6SHong Zhang /*@
1625814e85d6SHong Zhang  TSAdjointCostIntegral - Evaluate the cost integral in the adjoint run.
1626814e85d6SHong Zhang 
1627814e85d6SHong Zhang  Collective on TS
1628814e85d6SHong Zhang 
1629814e85d6SHong Zhang  Input Arguments:
1630814e85d6SHong Zhang  .  ts - time stepping context
1631814e85d6SHong Zhang 
1632814e85d6SHong Zhang  Level: advanced
1633814e85d6SHong Zhang 
1634814e85d6SHong Zhang  Notes:
1635814e85d6SHong Zhang  This function cannot be called until TSAdjointStep() has been completed.
1636814e85d6SHong Zhang 
1637814e85d6SHong Zhang  .seealso: TSAdjointSolve(), TSAdjointStep
1638814e85d6SHong Zhang  @*/
1639814e85d6SHong Zhang PetscErrorCode TSAdjointCostIntegral(TS ts)
1640814e85d6SHong Zhang {
1641814e85d6SHong Zhang     PetscErrorCode ierr;
1642814e85d6SHong Zhang     PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1643814e85d6SHong 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);
1644814e85d6SHong Zhang     ierr = (*ts->ops->adjointintegral)(ts);CHKERRQ(ierr);
1645814e85d6SHong Zhang     PetscFunctionReturn(0);
1646814e85d6SHong Zhang }
1647814e85d6SHong Zhang 
1648814e85d6SHong Zhang /* ------------------ Forward (tangent linear) sensitivity  ------------------*/
1649814e85d6SHong Zhang 
1650814e85d6SHong Zhang /*@
1651814e85d6SHong Zhang   TSForwardSetUp - Sets up the internal data structures for the later use
1652814e85d6SHong Zhang   of forward sensitivity analysis
1653814e85d6SHong Zhang 
1654814e85d6SHong Zhang   Collective on TS
1655814e85d6SHong Zhang 
1656814e85d6SHong Zhang   Input Parameter:
1657814e85d6SHong Zhang . ts - the TS context obtained from TSCreate()
1658814e85d6SHong Zhang 
1659814e85d6SHong Zhang   Level: advanced
1660814e85d6SHong Zhang 
1661814e85d6SHong Zhang .seealso: TSCreate(), TSDestroy(), TSSetUp()
1662814e85d6SHong Zhang @*/
1663814e85d6SHong Zhang PetscErrorCode TSForwardSetUp(TS ts)
1664814e85d6SHong Zhang {
1665814e85d6SHong Zhang   PetscErrorCode ierr;
1666814e85d6SHong Zhang 
1667814e85d6SHong Zhang   PetscFunctionBegin;
1668814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1669814e85d6SHong Zhang   if (ts->forwardsetupcalled) PetscFunctionReturn(0);
1670814e85d6SHong Zhang   if (ts->ops->forwardsetup) {
1671814e85d6SHong Zhang     ierr = (*ts->ops->forwardsetup)(ts);CHKERRQ(ierr);
1672814e85d6SHong Zhang   }
16737207e0fdSHong Zhang   ierr = VecDuplicate(ts->vec_sol,&ts->vec_sensip_col);CHKERRQ(ierr);
1674814e85d6SHong Zhang   ts->forwardsetupcalled = PETSC_TRUE;
1675814e85d6SHong Zhang   PetscFunctionReturn(0);
1676814e85d6SHong Zhang }
1677814e85d6SHong Zhang 
1678814e85d6SHong Zhang /*@
16797adebddeSHong Zhang   TSForwardReset - Reset the internal data structures used by forward sensitivity analysis
16807adebddeSHong Zhang 
16817adebddeSHong Zhang   Collective on TS
16827adebddeSHong Zhang 
16837adebddeSHong Zhang   Input Parameter:
16847adebddeSHong Zhang . ts - the TS context obtained from TSCreate()
16857adebddeSHong Zhang 
16867adebddeSHong Zhang   Level: advanced
16877adebddeSHong Zhang 
16887adebddeSHong Zhang .seealso: TSCreate(), TSDestroy(), TSForwardSetUp()
16897adebddeSHong Zhang @*/
16907adebddeSHong Zhang PetscErrorCode TSForwardReset(TS ts)
16917adebddeSHong Zhang {
16927207e0fdSHong Zhang   TS             quadts = ts->quadraturets;
16937adebddeSHong Zhang   PetscErrorCode ierr;
16947adebddeSHong Zhang 
16957adebddeSHong Zhang   PetscFunctionBegin;
16967adebddeSHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
16977adebddeSHong Zhang   if (ts->ops->forwardreset) {
16987adebddeSHong Zhang     ierr = (*ts->ops->forwardreset)(ts);CHKERRQ(ierr);
16997adebddeSHong Zhang   }
17007adebddeSHong Zhang   ierr = MatDestroy(&ts->mat_sensip);CHKERRQ(ierr);
17017207e0fdSHong Zhang   if (quadts) {
17027207e0fdSHong Zhang     ierr = MatDestroy(&quadts->mat_sensip);CHKERRQ(ierr);
17037207e0fdSHong Zhang   }
17047207e0fdSHong Zhang   ierr = VecDestroy(&ts->vec_sensip_col);CHKERRQ(ierr);
17057adebddeSHong Zhang   ts->forward_solve      = PETSC_FALSE;
17067adebddeSHong Zhang   ts->forwardsetupcalled = PETSC_FALSE;
17077adebddeSHong Zhang   PetscFunctionReturn(0);
17087adebddeSHong Zhang }
17097adebddeSHong Zhang 
17107adebddeSHong Zhang /*@
1711814e85d6SHong Zhang   TSForwardSetIntegralGradients - Set the vectors holding forward sensitivities of the integral term.
1712814e85d6SHong Zhang 
1713814e85d6SHong Zhang   Input Parameter:
1714a2b725a8SWilliam Gropp + ts- the TS context obtained from TSCreate()
1715814e85d6SHong Zhang . numfwdint- number of integrals
1716a2b725a8SWilliam Gropp - vp = the vectors containing the gradients for each integral w.r.t. parameters
1717814e85d6SHong Zhang 
17187207e0fdSHong Zhang   Level: deprecated
1719814e85d6SHong Zhang 
1720814e85d6SHong Zhang .seealso: TSForwardGetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep()
1721814e85d6SHong Zhang @*/
1722814e85d6SHong Zhang PetscErrorCode TSForwardSetIntegralGradients(TS ts,PetscInt numfwdint,Vec *vp)
1723814e85d6SHong Zhang {
1724814e85d6SHong Zhang   PetscFunctionBegin;
1725814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1726814e85d6SHong 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()");
1727814e85d6SHong Zhang   if (!ts->numcost) ts->numcost = numfwdint;
1728814e85d6SHong Zhang 
1729814e85d6SHong Zhang   ts->vecs_integral_sensip = vp;
1730814e85d6SHong Zhang   PetscFunctionReturn(0);
1731814e85d6SHong Zhang }
1732814e85d6SHong Zhang 
1733814e85d6SHong Zhang /*@
1734814e85d6SHong Zhang   TSForwardGetIntegralGradients - Returns the forward sensitivities ofthe integral term.
1735814e85d6SHong Zhang 
1736814e85d6SHong Zhang   Input Parameter:
1737814e85d6SHong Zhang . ts- the TS context obtained from TSCreate()
1738814e85d6SHong Zhang 
1739814e85d6SHong Zhang   Output Parameter:
1740814e85d6SHong Zhang . vp = the vectors containing the gradients for each integral w.r.t. parameters
1741814e85d6SHong Zhang 
17427207e0fdSHong Zhang   Level: deprecated
1743814e85d6SHong Zhang 
1744814e85d6SHong Zhang .seealso: TSForwardSetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep()
1745814e85d6SHong Zhang @*/
1746814e85d6SHong Zhang PetscErrorCode TSForwardGetIntegralGradients(TS ts,PetscInt *numfwdint,Vec **vp)
1747814e85d6SHong Zhang {
1748814e85d6SHong Zhang   PetscFunctionBegin;
1749814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1750814e85d6SHong Zhang   PetscValidPointer(vp,3);
1751814e85d6SHong Zhang   if (numfwdint) *numfwdint = ts->numcost;
1752814e85d6SHong Zhang   if (vp) *vp = ts->vecs_integral_sensip;
1753814e85d6SHong Zhang   PetscFunctionReturn(0);
1754814e85d6SHong Zhang }
1755814e85d6SHong Zhang 
1756814e85d6SHong Zhang /*@
1757814e85d6SHong Zhang   TSForwardStep - Compute the forward sensitivity for one time step.
1758814e85d6SHong Zhang 
1759814e85d6SHong Zhang   Collective on TS
1760814e85d6SHong Zhang 
1761814e85d6SHong Zhang   Input Arguments:
1762814e85d6SHong Zhang . ts - time stepping context
1763814e85d6SHong Zhang 
1764814e85d6SHong Zhang   Level: advanced
1765814e85d6SHong Zhang 
1766814e85d6SHong Zhang   Notes:
1767814e85d6SHong Zhang   This function cannot be called until TSStep() has been completed.
1768814e85d6SHong Zhang 
1769814e85d6SHong Zhang .seealso: TSForwardSetSensitivities(), TSForwardGetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardSetUp()
1770814e85d6SHong Zhang @*/
1771814e85d6SHong Zhang PetscErrorCode TSForwardStep(TS ts)
1772814e85d6SHong Zhang {
1773814e85d6SHong Zhang   PetscErrorCode ierr;
1774814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1775814e85d6SHong Zhang   if (!ts->ops->forwardstep) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"%s does not provide forward sensitivity analysis",((PetscObject)ts)->type_name);
1776814e85d6SHong Zhang   ierr = PetscLogEventBegin(TS_ForwardStep,ts,0,0,0);CHKERRQ(ierr);
1777814e85d6SHong Zhang   ierr = (*ts->ops->forwardstep)(ts);CHKERRQ(ierr);
1778814e85d6SHong Zhang   ierr = PetscLogEventEnd(TS_ForwardStep,ts,0,0,0);CHKERRQ(ierr);
1779814e85d6SHong Zhang   PetscFunctionReturn(0);
1780814e85d6SHong Zhang }
1781814e85d6SHong Zhang 
1782814e85d6SHong Zhang /*@
1783814e85d6SHong Zhang   TSForwardSetSensitivities - Sets the initial value of the trajectory sensitivities of solution  w.r.t. the problem parameters and initial values.
1784814e85d6SHong Zhang 
1785d083f849SBarry Smith   Logically Collective on TS
1786814e85d6SHong Zhang 
1787814e85d6SHong Zhang   Input Parameters:
1788814e85d6SHong Zhang + ts - the TS context obtained from TSCreate()
1789814e85d6SHong Zhang . nump - number of parameters
1790814e85d6SHong Zhang - Smat - sensitivities with respect to the parameters, the number of entries in these vectors is the same as the number of parameters
1791814e85d6SHong Zhang 
1792814e85d6SHong Zhang   Level: beginner
1793814e85d6SHong Zhang 
1794814e85d6SHong Zhang   Notes:
1795814e85d6SHong Zhang   Forward sensitivity is also called 'trajectory sensitivity' in some fields such as power systems.
1796814e85d6SHong Zhang   This function turns on a flag to trigger TSSolve() to compute forward sensitivities automatically.
1797814e85d6SHong Zhang   You must call this function before TSSolve().
1798814e85d6SHong Zhang   The entries in the sensitivity matrix must be correctly initialized with the values S = dy/dp|startingtime.
1799814e85d6SHong Zhang 
1800814e85d6SHong Zhang .seealso: TSForwardGetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep()
1801814e85d6SHong Zhang @*/
1802814e85d6SHong Zhang PetscErrorCode TSForwardSetSensitivities(TS ts,PetscInt nump,Mat Smat)
1803814e85d6SHong Zhang {
1804814e85d6SHong Zhang   PetscErrorCode ierr;
1805814e85d6SHong Zhang 
1806814e85d6SHong Zhang   PetscFunctionBegin;
1807814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1808814e85d6SHong Zhang   PetscValidHeaderSpecific(Smat,MAT_CLASSID,3);
1809814e85d6SHong Zhang   ts->forward_solve  = PETSC_TRUE;
1810b5b4017aSHong Zhang   if (nump == PETSC_DEFAULT) {
1811b5b4017aSHong Zhang     ierr = MatGetSize(Smat,NULL,&ts->num_parameters);CHKERRQ(ierr);
1812b5b4017aSHong Zhang   } else ts->num_parameters = nump;
1813814e85d6SHong Zhang   ierr = PetscObjectReference((PetscObject)Smat);CHKERRQ(ierr);
1814814e85d6SHong Zhang   ierr = MatDestroy(&ts->mat_sensip);CHKERRQ(ierr);
1815814e85d6SHong Zhang   ts->mat_sensip = Smat;
1816814e85d6SHong Zhang   PetscFunctionReturn(0);
1817814e85d6SHong Zhang }
1818814e85d6SHong Zhang 
1819814e85d6SHong Zhang /*@
1820814e85d6SHong Zhang   TSForwardGetSensitivities - Returns the trajectory sensitivities
1821814e85d6SHong Zhang 
1822814e85d6SHong Zhang   Not Collective, but Vec returned is parallel if TS is parallel
1823814e85d6SHong Zhang 
1824814e85d6SHong Zhang   Output Parameter:
1825814e85d6SHong Zhang + ts - the TS context obtained from TSCreate()
1826814e85d6SHong Zhang . nump - number of parameters
1827814e85d6SHong Zhang - Smat - sensitivities with respect to the parameters, the number of entries in these vectors is the same as the number of parameters
1828814e85d6SHong Zhang 
1829814e85d6SHong Zhang   Level: intermediate
1830814e85d6SHong Zhang 
1831814e85d6SHong Zhang .seealso: TSForwardSetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep()
1832814e85d6SHong Zhang @*/
1833814e85d6SHong Zhang PetscErrorCode TSForwardGetSensitivities(TS ts,PetscInt *nump,Mat *Smat)
1834814e85d6SHong Zhang {
1835814e85d6SHong Zhang   PetscFunctionBegin;
1836814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1837814e85d6SHong Zhang   if (nump) *nump = ts->num_parameters;
1838814e85d6SHong Zhang   if (Smat) *Smat = ts->mat_sensip;
1839814e85d6SHong Zhang   PetscFunctionReturn(0);
1840814e85d6SHong Zhang }
1841814e85d6SHong Zhang 
1842814e85d6SHong Zhang /*@
1843814e85d6SHong Zhang    TSForwardCostIntegral - Evaluate the cost integral in the forward run.
1844814e85d6SHong Zhang 
1845814e85d6SHong Zhang    Collective on TS
1846814e85d6SHong Zhang 
1847814e85d6SHong Zhang    Input Arguments:
1848814e85d6SHong Zhang .  ts - time stepping context
1849814e85d6SHong Zhang 
1850814e85d6SHong Zhang    Level: advanced
1851814e85d6SHong Zhang 
1852814e85d6SHong Zhang    Notes:
1853814e85d6SHong Zhang    This function cannot be called until TSStep() has been completed.
1854814e85d6SHong Zhang 
1855814e85d6SHong Zhang .seealso: TSSolve(), TSAdjointCostIntegral()
1856814e85d6SHong Zhang @*/
1857814e85d6SHong Zhang PetscErrorCode TSForwardCostIntegral(TS ts)
1858814e85d6SHong Zhang {
1859814e85d6SHong Zhang   PetscErrorCode ierr;
1860814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1861814e85d6SHong 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);
1862814e85d6SHong Zhang   ierr = (*ts->ops->forwardintegral)(ts);CHKERRQ(ierr);
1863814e85d6SHong Zhang   PetscFunctionReturn(0);
1864814e85d6SHong Zhang }
1865b5b4017aSHong Zhang 
1866b5b4017aSHong Zhang /*@
1867b5b4017aSHong Zhang   TSForwardSetInitialSensitivities - Set initial values for tangent linear sensitivities
1868b5b4017aSHong Zhang 
1869d083f849SBarry Smith   Collective on TS
1870b5b4017aSHong Zhang 
1871b5b4017aSHong Zhang   Input Parameter
1872b5b4017aSHong Zhang + ts - the TS context obtained from TSCreate()
1873b5b4017aSHong Zhang - didp - parametric sensitivities of the initial condition
1874b5b4017aSHong Zhang 
1875b5b4017aSHong Zhang   Level: intermediate
1876b5b4017aSHong Zhang 
1877b5b4017aSHong 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.
1878b5b4017aSHong Zhang 
1879b5b4017aSHong Zhang .seealso: TSForwardSetSensitivities()
1880b5b4017aSHong Zhang @*/
1881b5b4017aSHong Zhang PetscErrorCode TSForwardSetInitialSensitivities(TS ts,Mat didp)
1882b5b4017aSHong Zhang {
1883b5b4017aSHong Zhang   PetscErrorCode ierr;
1884b5b4017aSHong Zhang 
1885b5b4017aSHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1886b5b4017aSHong Zhang   PetscValidHeaderSpecific(didp,MAT_CLASSID,2);
1887b5b4017aSHong Zhang   if (!ts->mat_sensip) {
1888b5b4017aSHong Zhang     ierr = TSForwardSetSensitivities(ts,PETSC_DEFAULT,didp);CHKERRQ(ierr);
1889b5b4017aSHong Zhang   }
1890b5b4017aSHong Zhang   PetscFunctionReturn(0);
1891b5b4017aSHong Zhang }
18926affb6f8SHong Zhang 
18936affb6f8SHong Zhang /*@
18946affb6f8SHong Zhang    TSForwardGetStages - Get the number of stages and the tangent linear sensitivities at the intermediate stages
18956affb6f8SHong Zhang 
18966affb6f8SHong Zhang    Input Parameter:
18976affb6f8SHong Zhang .  ts - the TS context obtained from TSCreate()
18986affb6f8SHong Zhang 
18996affb6f8SHong Zhang    Output Parameters:
1900cd4cee2dSHong Zhang +  ns - number of stages
19016affb6f8SHong Zhang -  S - tangent linear sensitivities at the intermediate stages
19026affb6f8SHong Zhang 
19036affb6f8SHong Zhang    Level: advanced
19046affb6f8SHong Zhang 
19056affb6f8SHong Zhang @*/
19066affb6f8SHong Zhang PetscErrorCode TSForwardGetStages(TS ts,PetscInt *ns,Mat **S)
19076affb6f8SHong Zhang {
19086affb6f8SHong Zhang   PetscErrorCode ierr;
19096affb6f8SHong Zhang 
19106affb6f8SHong Zhang   PetscFunctionBegin;
19116affb6f8SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
19126affb6f8SHong Zhang 
19136affb6f8SHong Zhang   if (!ts->ops->getstages) *S=NULL;
19146affb6f8SHong Zhang   else {
19156affb6f8SHong Zhang     ierr = (*ts->ops->forwardgetstages)(ts,ns,S);CHKERRQ(ierr);
19166affb6f8SHong Zhang   }
19176affb6f8SHong Zhang   PetscFunctionReturn(0);
19186affb6f8SHong Zhang }
1919cd4cee2dSHong Zhang 
1920cd4cee2dSHong Zhang /*@
1921cd4cee2dSHong Zhang    TSCreateQuadratureTS - Create a sub-TS that evaluates integrals over time
1922cd4cee2dSHong Zhang 
1923cd4cee2dSHong Zhang    Input Parameter:
1924cd4cee2dSHong Zhang +  ts - the TS context obtained from TSCreate()
1925cd4cee2dSHong Zhang -  fwd - flag indicating whether to evaluate cost integral in the forward run or the adjoint run
1926cd4cee2dSHong Zhang 
1927cd4cee2dSHong Zhang    Output Parameters:
1928cd4cee2dSHong Zhang .  quadts - the child TS context
1929cd4cee2dSHong Zhang 
1930cd4cee2dSHong Zhang    Level: intermediate
1931cd4cee2dSHong Zhang 
1932cd4cee2dSHong Zhang .seealso: TSGetQuadratureTS()
1933cd4cee2dSHong Zhang @*/
1934cd4cee2dSHong Zhang PetscErrorCode TSCreateQuadratureTS(TS ts,PetscBool fwd,TS *quadts)
1935cd4cee2dSHong Zhang {
1936cd4cee2dSHong Zhang   char prefix[128];
1937cd4cee2dSHong Zhang   PetscErrorCode ierr;
1938cd4cee2dSHong Zhang 
1939cd4cee2dSHong Zhang   PetscFunctionBegin;
1940cd4cee2dSHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1941cd4cee2dSHong Zhang   PetscValidPointer(quadts,2);
1942cd4cee2dSHong Zhang   ierr = TSDestroy(&ts->quadraturets);CHKERRQ(ierr);
1943cd4cee2dSHong Zhang   ierr = TSCreate(PetscObjectComm((PetscObject)ts),&ts->quadraturets);CHKERRQ(ierr);
1944cd4cee2dSHong Zhang   ierr = PetscObjectIncrementTabLevel((PetscObject)ts->quadraturets,(PetscObject)ts,1);CHKERRQ(ierr);
1945cd4cee2dSHong Zhang   ierr = PetscLogObjectParent((PetscObject)ts,(PetscObject)ts->quadraturets);CHKERRQ(ierr);
1946cd4cee2dSHong Zhang   ierr = PetscSNPrintf(prefix,sizeof(prefix),"%squad_",((PetscObject)ts)->prefix ? ((PetscObject)ts)->prefix : "");CHKERRQ(ierr);
1947cd4cee2dSHong Zhang   ierr = TSSetOptionsPrefix(ts->quadraturets,prefix);CHKERRQ(ierr);
1948cd4cee2dSHong Zhang   *quadts = ts->quadraturets;
1949cd4cee2dSHong Zhang 
1950cd4cee2dSHong Zhang   if (ts->numcost) {
1951cd4cee2dSHong Zhang     ierr = VecCreateSeq(PETSC_COMM_SELF,ts->numcost,&(*quadts)->vec_sol);CHKERRQ(ierr);
1952cd4cee2dSHong Zhang   } else {
1953cd4cee2dSHong Zhang     ierr = VecCreateSeq(PETSC_COMM_SELF,1,&(*quadts)->vec_sol);CHKERRQ(ierr);
1954cd4cee2dSHong Zhang   }
1955cd4cee2dSHong Zhang   ts->costintegralfwd = fwd;
1956cd4cee2dSHong Zhang   PetscFunctionReturn(0);
1957cd4cee2dSHong Zhang }
1958cd4cee2dSHong Zhang 
1959cd4cee2dSHong Zhang /*@
1960cd4cee2dSHong Zhang    TSGetQuadratureTS - Return the sub-TS that evaluates integrals over time
1961cd4cee2dSHong Zhang 
1962cd4cee2dSHong Zhang    Input Parameter:
1963cd4cee2dSHong Zhang .  ts - the TS context obtained from TSCreate()
1964cd4cee2dSHong Zhang 
1965cd4cee2dSHong Zhang    Output Parameters:
1966cd4cee2dSHong Zhang +  fwd - flag indicating whether to evaluate cost integral in the forward run or the adjoint run
1967cd4cee2dSHong Zhang -  quadts - the child TS context
1968cd4cee2dSHong Zhang 
1969cd4cee2dSHong Zhang    Level: intermediate
1970cd4cee2dSHong Zhang 
1971cd4cee2dSHong Zhang .seealso: TSCreateQuadratureTS()
1972cd4cee2dSHong Zhang @*/
1973cd4cee2dSHong Zhang PetscErrorCode TSGetQuadratureTS(TS ts,PetscBool *fwd,TS *quadts)
1974cd4cee2dSHong Zhang {
1975cd4cee2dSHong Zhang   PetscFunctionBegin;
1976cd4cee2dSHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1977cd4cee2dSHong Zhang   if (fwd) *fwd = ts->costintegralfwd;
1978cd4cee2dSHong Zhang   if (quadts) *quadts = ts->quadraturets;
1979cd4cee2dSHong Zhang   PetscFunctionReturn(0);
1980cd4cee2dSHong Zhang }
1981*ba0675f6SHong Zhang 
1982*ba0675f6SHong Zhang /*@
1983*ba0675f6SHong Zhang    TSComputeSNESJacobian - Compute the SNESJacobian
1984*ba0675f6SHong Zhang 
1985*ba0675f6SHong Zhang    Input Parameters:
1986*ba0675f6SHong Zhang +  ts - the TS context obtained from TSCreate()
1987*ba0675f6SHong Zhang -  x - state vector
1988*ba0675f6SHong Zhang 
1989*ba0675f6SHong Zhang    Output Parameters:
1990*ba0675f6SHong Zhang +  J - Jacobian matrix
1991*ba0675f6SHong Zhang -  Jpre - preconditioning matrix for J (may be same as J)
1992*ba0675f6SHong Zhang 
1993*ba0675f6SHong Zhang    Level: developer
1994*ba0675f6SHong Zhang 
1995*ba0675f6SHong Zhang    Notes:
1996*ba0675f6SHong Zhang    Using SNES to compute the Jacobian enables finite differencing when TS Jacobian is not available.
1997*ba0675f6SHong Zhang @*/
1998*ba0675f6SHong Zhang PetscErrorCode TSComputeSNESJacobian(TS ts,Vec x,Mat J,Mat Jpre)
1999*ba0675f6SHong Zhang {
2000*ba0675f6SHong Zhang   SNES           snes = ts->snes;
2001*ba0675f6SHong Zhang   PetscErrorCode (*jac)(SNES,Vec,Mat,Mat,void*) = NULL;
2002*ba0675f6SHong Zhang   PetscErrorCode ierr;
2003*ba0675f6SHong Zhang 
2004*ba0675f6SHong Zhang   PetscFunctionBegin;
2005*ba0675f6SHong Zhang   /*
2006*ba0675f6SHong Zhang     Unlike implicit methods, explicit methods do not have SNESMatFDColoring in the snes object
2007*ba0675f6SHong Zhang     because SNESSolve() has not been called yet; so querying SNESMatFDColoring does not work for
2008*ba0675f6SHong Zhang     explicit methods. Instead, we check the Jacobian compute function directly to determin if FD
2009*ba0675f6SHong Zhang     coloring is used.
2010*ba0675f6SHong Zhang   */
2011*ba0675f6SHong Zhang   ierr = SNESGetJacobian(snes,NULL,NULL,&jac,NULL);CHKERRQ(ierr);
2012*ba0675f6SHong Zhang   if (jac == SNESComputeJacobianDefaultColor) {
2013*ba0675f6SHong Zhang     Vec f;
2014*ba0675f6SHong Zhang     ierr = SNESSetSolution(snes,x);CHKERRQ(ierr);
2015*ba0675f6SHong Zhang     ierr = SNESGetFunction(snes,&f,NULL,NULL);CHKERRQ(ierr);
2016*ba0675f6SHong Zhang     /* Force MatFDColoringApply to evaluate the SNES residual function for the base vector */
2017*ba0675f6SHong Zhang     ierr = SNESComputeFunction(snes,x,f);CHKERRQ(ierr);
2018*ba0675f6SHong Zhang   }
2019*ba0675f6SHong Zhang   ierr = SNESComputeJacobian(snes,x,J,Jpre);CHKERRQ(ierr);
2020*ba0675f6SHong Zhang   PetscFunctionReturn(0);
2021*ba0675f6SHong Zhang }
2022