xref: /petsc/src/ts/interface/sensitivity/tssen.c (revision 80ab5e31fb2bb9538e7952a0ebdda3c9178394af)
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 
13c3339decSBarry Smith   Logically Collective
14814e85d6SHong Zhang 
15814e85d6SHong Zhang   Input Parameters:
16bcf0153eSBarry Smith + 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 
30bcf0153eSBarry Smith   Note:
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 
33bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSGetRHSJacobianP()`
34814e85d6SHong Zhang @*/
35d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSetRHSJacobianP(TS ts, Mat Amat, PetscErrorCode (*func)(TS, PetscReal, Vec, Mat, void *), void *ctx)
36d71ae5a4SJacob Faibussowitsch {
37814e85d6SHong Zhang   PetscFunctionBegin;
38814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
39814e85d6SHong Zhang   PetscValidHeaderSpecific(Amat, MAT_CLASSID, 2);
40814e85d6SHong Zhang 
41814e85d6SHong Zhang   ts->rhsjacobianp    = func;
42814e85d6SHong Zhang   ts->rhsjacobianpctx = ctx;
43814e85d6SHong Zhang   if (Amat) {
449566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)Amat));
459566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&ts->Jacprhs));
4633f72c5dSHong Zhang     ts->Jacprhs = Amat;
47814e85d6SHong Zhang   }
483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
49814e85d6SHong Zhang }
50814e85d6SHong Zhang 
51814e85d6SHong Zhang /*@C
52cd4cee2dSHong 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.
53cd4cee2dSHong Zhang 
54c3339decSBarry Smith   Logically Collective
55cd4cee2dSHong Zhang 
56f899ff85SJose E. Roman   Input Parameter:
57bcf0153eSBarry Smith . ts - `TS` context obtained from `TSCreate()`
58cd4cee2dSHong Zhang 
59cd4cee2dSHong Zhang   Output Parameters:
60cd4cee2dSHong Zhang + Amat - JacobianP matrix
61cd4cee2dSHong Zhang . func - function
62cd4cee2dSHong Zhang - ctx - [optional] user-defined function context
63cd4cee2dSHong Zhang 
64cd4cee2dSHong Zhang   Calling sequence of func:
65cd4cee2dSHong Zhang $ func (TS ts,PetscReal t,Vec y,Mat A,void *ctx);
66cd4cee2dSHong Zhang +   t - current timestep
67cd4cee2dSHong Zhang .   U - input vector (current ODE solution)
68cd4cee2dSHong Zhang .   A - output matrix
69cd4cee2dSHong Zhang -   ctx - [optional] user-defined function context
70cd4cee2dSHong Zhang 
71cd4cee2dSHong Zhang   Level: intermediate
72cd4cee2dSHong Zhang 
73bcf0153eSBarry Smith   Note:
74cd4cee2dSHong 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
75cd4cee2dSHong Zhang 
76bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSSetRHSJacobianP()`, `TS`, `TSGetRHSJacobianP()`
77cd4cee2dSHong Zhang @*/
78d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGetRHSJacobianP(TS ts, Mat *Amat, PetscErrorCode (**func)(TS, PetscReal, Vec, Mat, void *), void **ctx)
79d71ae5a4SJacob Faibussowitsch {
80cd4cee2dSHong Zhang   PetscFunctionBegin;
81cd4cee2dSHong Zhang   if (func) *func = ts->rhsjacobianp;
82cd4cee2dSHong Zhang   if (ctx) *ctx = ts->rhsjacobianpctx;
83cd4cee2dSHong Zhang   if (Amat) *Amat = ts->Jacprhs;
843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
85cd4cee2dSHong Zhang }
86cd4cee2dSHong Zhang 
87cd4cee2dSHong Zhang /*@C
88814e85d6SHong Zhang   TSComputeRHSJacobianP - Runs the user-defined JacobianP function.
89814e85d6SHong Zhang 
90c3339decSBarry Smith   Collective
91814e85d6SHong Zhang 
92814e85d6SHong Zhang   Input Parameters:
93bcf0153eSBarry Smith . ts   - The `TS` context obtained from `TSCreate()`
94814e85d6SHong Zhang 
95814e85d6SHong Zhang   Level: developer
96814e85d6SHong Zhang 
97bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSSetRHSJacobianP()`, `TS`
98814e85d6SHong Zhang @*/
99d71ae5a4SJacob Faibussowitsch PetscErrorCode TSComputeRHSJacobianP(TS ts, PetscReal t, Vec U, Mat Amat)
100d71ae5a4SJacob Faibussowitsch {
101814e85d6SHong Zhang   PetscFunctionBegin;
1023ba16761SJacob Faibussowitsch   if (!Amat) PetscFunctionReturn(PETSC_SUCCESS);
103814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
104c9ad9de0SHong Zhang   PetscValidHeaderSpecific(U, VEC_CLASSID, 3);
105814e85d6SHong Zhang 
106*80ab5e31SHong Zhang   if (ts->rhsjacobianp) PetscCallBack("TS callback JacobianP for sensitivity analysis", (*ts->rhsjacobianp)(ts, t, U, Amat, ts->rhsjacobianpctx));
107*80ab5e31SHong Zhang   else {
108*80ab5e31SHong Zhang     PetscBool assembled;
109*80ab5e31SHong Zhang     PetscCall(MatZeroEntries(Amat));
110*80ab5e31SHong Zhang     PetscCall(MatAssembled(Amat, &assembled));
111*80ab5e31SHong Zhang     if (!assembled) {
112*80ab5e31SHong Zhang       PetscCall(MatAssemblyBegin(Amat, MAT_FINAL_ASSEMBLY));
113*80ab5e31SHong Zhang       PetscCall(MatAssemblyEnd(Amat, MAT_FINAL_ASSEMBLY));
114*80ab5e31SHong Zhang     }
115*80ab5e31SHong Zhang   }
1163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
117814e85d6SHong Zhang }
118814e85d6SHong Zhang 
119814e85d6SHong Zhang /*@C
12033f72c5dSHong 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.
12133f72c5dSHong Zhang 
122c3339decSBarry Smith   Logically Collective
12333f72c5dSHong Zhang 
12433f72c5dSHong Zhang   Input Parameters:
125bcf0153eSBarry Smith + ts - `TS` context obtained from `TSCreate()`
12633f72c5dSHong Zhang . Amat - JacobianP matrix
12733f72c5dSHong Zhang . func - function
12833f72c5dSHong Zhang - ctx - [optional] user-defined function context
12933f72c5dSHong Zhang 
13033f72c5dSHong Zhang   Calling sequence of func:
13133f72c5dSHong Zhang $ func (TS ts,PetscReal t,Vec y,Mat A,void *ctx);
13233f72c5dSHong Zhang +   t - current timestep
13333f72c5dSHong Zhang .   U - input vector (current ODE solution)
13433f72c5dSHong Zhang .   Udot - time derivative of state vector
13533f72c5dSHong Zhang .   shift - shift to apply, see note below
13633f72c5dSHong Zhang .   A - output matrix
13733f72c5dSHong Zhang -   ctx - [optional] user-defined function context
13833f72c5dSHong Zhang 
13933f72c5dSHong Zhang   Level: intermediate
14033f72c5dSHong Zhang 
141bcf0153eSBarry Smith   Note:
14233f72c5dSHong 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
14333f72c5dSHong Zhang 
144bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSSetRHSJacobianP()`, `TS`
14533f72c5dSHong Zhang @*/
146d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSetIJacobianP(TS ts, Mat Amat, PetscErrorCode (*func)(TS, PetscReal, Vec, Vec, PetscReal, Mat, void *), void *ctx)
147d71ae5a4SJacob Faibussowitsch {
14833f72c5dSHong Zhang   PetscFunctionBegin;
14933f72c5dSHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
15033f72c5dSHong Zhang   PetscValidHeaderSpecific(Amat, MAT_CLASSID, 2);
15133f72c5dSHong Zhang 
15233f72c5dSHong Zhang   ts->ijacobianp    = func;
15333f72c5dSHong Zhang   ts->ijacobianpctx = ctx;
15433f72c5dSHong Zhang   if (Amat) {
1559566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)Amat));
1569566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&ts->Jacp));
15733f72c5dSHong Zhang     ts->Jacp = Amat;
15833f72c5dSHong Zhang   }
1593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16033f72c5dSHong Zhang }
16133f72c5dSHong Zhang 
16233f72c5dSHong Zhang /*@C
16333f72c5dSHong Zhang   TSComputeIJacobianP - Runs the user-defined IJacobianP function.
16433f72c5dSHong Zhang 
165c3339decSBarry Smith   Collective
16633f72c5dSHong Zhang 
16733f72c5dSHong Zhang   Input Parameters:
168bcf0153eSBarry Smith + ts - the `TS` context
16933f72c5dSHong Zhang . t - current timestep
17033f72c5dSHong Zhang . U - state vector
17133f72c5dSHong Zhang . Udot - time derivative of state vector
17233f72c5dSHong Zhang . shift - shift to apply, see note below
173*80ab5e31SHong Zhang - imex - flag indicates if the method is IMEX so that the RHSJacobianP should be kept separate
17433f72c5dSHong Zhang 
17533f72c5dSHong Zhang   Output Parameters:
17633f72c5dSHong Zhang . A - Jacobian matrix
17733f72c5dSHong Zhang 
17833f72c5dSHong Zhang   Level: developer
17933f72c5dSHong Zhang 
180bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSSetIJacobianP()`
18133f72c5dSHong Zhang @*/
182d71ae5a4SJacob Faibussowitsch PetscErrorCode TSComputeIJacobianP(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal shift, Mat Amat, PetscBool imex)
183d71ae5a4SJacob Faibussowitsch {
18433f72c5dSHong Zhang   PetscFunctionBegin;
1853ba16761SJacob Faibussowitsch   if (!Amat) PetscFunctionReturn(PETSC_SUCCESS);
18633f72c5dSHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
18733f72c5dSHong Zhang   PetscValidHeaderSpecific(U, VEC_CLASSID, 3);
18833f72c5dSHong Zhang   PetscValidHeaderSpecific(Udot, VEC_CLASSID, 4);
18933f72c5dSHong Zhang 
1909566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(TS_JacobianPEval, ts, U, Amat, 0));
19148a46eb9SPierre Jolivet   if (ts->ijacobianp) PetscCallBack("TS callback JacobianP for sensitivity analysis", (*ts->ijacobianp)(ts, t, U, Udot, shift, Amat, ts->ijacobianpctx));
19233f72c5dSHong Zhang   if (imex) {
19333f72c5dSHong Zhang     if (!ts->ijacobianp) { /* system was written as Udot = G(t,U) */
19433f72c5dSHong Zhang       PetscBool assembled;
1959566063dSJacob Faibussowitsch       PetscCall(MatZeroEntries(Amat));
1969566063dSJacob Faibussowitsch       PetscCall(MatAssembled(Amat, &assembled));
19733f72c5dSHong Zhang       if (!assembled) {
1989566063dSJacob Faibussowitsch         PetscCall(MatAssemblyBegin(Amat, MAT_FINAL_ASSEMBLY));
1999566063dSJacob Faibussowitsch         PetscCall(MatAssemblyEnd(Amat, MAT_FINAL_ASSEMBLY));
20033f72c5dSHong Zhang       }
20133f72c5dSHong Zhang     }
20233f72c5dSHong Zhang   } else {
2031baa6e33SBarry Smith     if (ts->rhsjacobianp) PetscCall(TSComputeRHSJacobianP(ts, t, U, ts->Jacprhs));
20433f72c5dSHong Zhang     if (ts->Jacprhs == Amat) { /* No IJacobian, so we only have the RHS matrix */
2059566063dSJacob Faibussowitsch       PetscCall(MatScale(Amat, -1));
20633f72c5dSHong Zhang     } else if (ts->Jacprhs) { /* Both IJacobian and RHSJacobian */
20733f72c5dSHong Zhang       MatStructure axpy = DIFFERENT_NONZERO_PATTERN;
20833f72c5dSHong Zhang       if (!ts->ijacobianp) { /* No IJacobianp provided, but we have a separate RHS matrix */
2099566063dSJacob Faibussowitsch         PetscCall(MatZeroEntries(Amat));
21033f72c5dSHong Zhang       }
2119566063dSJacob Faibussowitsch       PetscCall(MatAXPY(Amat, -1, ts->Jacprhs, axpy));
21233f72c5dSHong Zhang     }
21333f72c5dSHong Zhang   }
2149566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(TS_JacobianPEval, ts, U, Amat, 0));
2153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
21633f72c5dSHong Zhang }
21733f72c5dSHong Zhang 
21833f72c5dSHong Zhang /*@C
219814e85d6SHong Zhang     TSSetCostIntegrand - Sets the routine for evaluating the integral term in one or more cost functions
220814e85d6SHong Zhang 
221c3339decSBarry Smith     Logically Collective
222814e85d6SHong Zhang 
223814e85d6SHong Zhang     Input Parameters:
224bcf0153eSBarry Smith +   ts - the `TS` context obtained from `TSCreate()`
225814e85d6SHong Zhang .   numcost - number of gradients to be computed, this is the number of cost functions
226814e85d6SHong Zhang .   costintegral - vector that stores the integral values
227814e85d6SHong Zhang .   rf - routine for evaluating the integrand function
228c9ad9de0SHong Zhang .   drduf - function that computes the gradients of the r's with respect to u
229814e85d6SHong 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)
230814e85d6SHong Zhang .   fwd - flag indicating whether to evaluate cost integral in the forward run or the adjoint run
231814e85d6SHong Zhang -   ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
232814e85d6SHong Zhang 
233814e85d6SHong Zhang     Calling sequence of rf:
234c9ad9de0SHong Zhang $   PetscErrorCode rf(TS ts,PetscReal t,Vec U,Vec F,void *ctx);
235814e85d6SHong Zhang 
236c9ad9de0SHong Zhang     Calling sequence of drduf:
237c9ad9de0SHong Zhang $   PetscErroCode drduf(TS ts,PetscReal t,Vec U,Vec *dRdU,void *ctx);
238814e85d6SHong Zhang 
239814e85d6SHong Zhang     Calling sequence of drdpf:
240c9ad9de0SHong Zhang $   PetscErroCode drdpf(TS ts,PetscReal t,Vec U,Vec *dRdP,void *ctx);
241814e85d6SHong Zhang 
242cd4cee2dSHong Zhang     Level: deprecated
243814e85d6SHong Zhang 
244bcf0153eSBarry Smith     Note:
245814e85d6SHong Zhang     For optimization there is usually a single cost function (numcost = 1). For sensitivities there may be multiple cost functions
246814e85d6SHong Zhang 
247bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSSetRHSJacobianP()`, `TSGetCostGradients()`, `TSSetCostGradients()`
248814e85d6SHong Zhang @*/
249d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSetCostIntegrand(TS ts, PetscInt numcost, Vec costintegral, PetscErrorCode (*rf)(TS, PetscReal, Vec, Vec, void *), PetscErrorCode (*drduf)(TS, PetscReal, Vec, Vec *, void *), PetscErrorCode (*drdpf)(TS, PetscReal, Vec, Vec *, void *), PetscBool fwd, void *ctx)
250d71ae5a4SJacob Faibussowitsch {
251814e85d6SHong Zhang   PetscFunctionBegin;
252814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
253814e85d6SHong Zhang   if (costintegral) PetscValidHeaderSpecific(costintegral, VEC_CLASSID, 3);
2543c633725SBarry Smith   PetscCheck(!ts->numcost || ts->numcost == numcost, PetscObjectComm((PetscObject)ts), PETSC_ERR_USER, "The number of cost functions (2nd parameter of TSSetCostIntegrand()) is inconsistent with the one set by TSSetCostGradients() or TSForwardSetIntegralGradients()");
255814e85d6SHong Zhang   if (!ts->numcost) ts->numcost = numcost;
256814e85d6SHong Zhang 
257814e85d6SHong Zhang   if (costintegral) {
2589566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)costintegral));
2599566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&ts->vec_costintegral));
260814e85d6SHong Zhang     ts->vec_costintegral = costintegral;
261814e85d6SHong Zhang   } else {
262814e85d6SHong Zhang     if (!ts->vec_costintegral) { /* Create a seq vec if user does not provide one */
2639566063dSJacob Faibussowitsch       PetscCall(VecCreateSeq(PETSC_COMM_SELF, numcost, &ts->vec_costintegral));
264814e85d6SHong Zhang     } else {
2659566063dSJacob Faibussowitsch       PetscCall(VecSet(ts->vec_costintegral, 0.0));
266814e85d6SHong Zhang     }
267814e85d6SHong Zhang   }
268814e85d6SHong Zhang   if (!ts->vec_costintegrand) {
2699566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(ts->vec_costintegral, &ts->vec_costintegrand));
270814e85d6SHong Zhang   } else {
2719566063dSJacob Faibussowitsch     PetscCall(VecSet(ts->vec_costintegrand, 0.0));
272814e85d6SHong Zhang   }
273814e85d6SHong Zhang   ts->costintegralfwd  = fwd; /* Evaluate the cost integral in forward run if fwd is true */
274814e85d6SHong Zhang   ts->costintegrand    = rf;
275814e85d6SHong Zhang   ts->costintegrandctx = ctx;
276c9ad9de0SHong Zhang   ts->drdufunction     = drduf;
277814e85d6SHong Zhang   ts->drdpfunction     = drdpf;
2783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
279814e85d6SHong Zhang }
280814e85d6SHong Zhang 
281b5b4017aSHong Zhang /*@C
282814e85d6SHong Zhang    TSGetCostIntegral - Returns the values of the integral term in the cost functions.
283814e85d6SHong Zhang    It is valid to call the routine after a backward run.
284814e85d6SHong Zhang 
285814e85d6SHong Zhang    Not Collective
286814e85d6SHong Zhang 
287814e85d6SHong Zhang    Input Parameter:
288bcf0153eSBarry Smith .  ts - the `TS` context obtained from `TSCreate()`
289814e85d6SHong Zhang 
290814e85d6SHong Zhang    Output Parameter:
291814e85d6SHong Zhang .  v - the vector containing the integrals for each cost function
292814e85d6SHong Zhang 
293814e85d6SHong Zhang    Level: intermediate
294814e85d6SHong Zhang 
295bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSAdjointSolve()`, ``TSSetCostIntegrand()`
296814e85d6SHong Zhang @*/
297d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGetCostIntegral(TS ts, Vec *v)
298d71ae5a4SJacob Faibussowitsch {
299cd4cee2dSHong Zhang   TS quadts;
300cd4cee2dSHong Zhang 
301814e85d6SHong Zhang   PetscFunctionBegin;
302814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
303814e85d6SHong Zhang   PetscValidPointer(v, 2);
3049566063dSJacob Faibussowitsch   PetscCall(TSGetQuadratureTS(ts, NULL, &quadts));
305cd4cee2dSHong Zhang   *v = quadts->vec_sol;
3063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
307814e85d6SHong Zhang }
308814e85d6SHong Zhang 
309b5b4017aSHong Zhang /*@C
310814e85d6SHong Zhang    TSComputeCostIntegrand - Evaluates the integral function in the cost functions.
311814e85d6SHong Zhang 
312814e85d6SHong Zhang    Input Parameters:
313bcf0153eSBarry Smith +  ts - the `TS` context
314814e85d6SHong Zhang .  t - current time
315c9ad9de0SHong Zhang -  U - state vector, i.e. current solution
316814e85d6SHong Zhang 
317814e85d6SHong Zhang    Output Parameter:
318c9ad9de0SHong Zhang .  Q - vector of size numcost to hold the outputs
319814e85d6SHong Zhang 
320bcf0153eSBarry Smith    Level: deprecated
321bcf0153eSBarry Smith 
322bcf0153eSBarry Smith    Note:
323814e85d6SHong Zhang    Most users should not need to explicitly call this routine, as it
324814e85d6SHong Zhang    is used internally within the sensitivity analysis context.
325814e85d6SHong Zhang 
326bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSAdjointSolve()`, `TSSetCostIntegrand()`
327814e85d6SHong Zhang @*/
328d71ae5a4SJacob Faibussowitsch PetscErrorCode TSComputeCostIntegrand(TS ts, PetscReal t, Vec U, Vec Q)
329d71ae5a4SJacob Faibussowitsch {
330814e85d6SHong Zhang   PetscFunctionBegin;
331814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
332c9ad9de0SHong Zhang   PetscValidHeaderSpecific(U, VEC_CLASSID, 3);
333c9ad9de0SHong Zhang   PetscValidHeaderSpecific(Q, VEC_CLASSID, 4);
334814e85d6SHong Zhang 
3359566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(TS_FunctionEval, ts, U, Q, 0));
336792fecdfSBarry Smith   if (ts->costintegrand) PetscCallBack("TS callback integrand in the cost function", (*ts->costintegrand)(ts, t, U, Q, ts->costintegrandctx));
337ef1023bdSBarry Smith   else PetscCall(VecZeroEntries(Q));
3389566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(TS_FunctionEval, ts, U, Q, 0));
3393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
340814e85d6SHong Zhang }
341814e85d6SHong Zhang 
342b5b4017aSHong Zhang /*@C
343bcf0153eSBarry Smith   TSComputeDRDUFunction - Deprecated, use `TSGetQuadratureTS()` then `TSComputeRHSJacobian()`
344814e85d6SHong Zhang 
345d76be551SHong Zhang   Level: deprecated
346814e85d6SHong Zhang 
347814e85d6SHong Zhang @*/
348d71ae5a4SJacob Faibussowitsch PetscErrorCode TSComputeDRDUFunction(TS ts, PetscReal t, Vec U, Vec *DRDU)
349d71ae5a4SJacob Faibussowitsch {
350814e85d6SHong Zhang   PetscFunctionBegin;
3513ba16761SJacob Faibussowitsch   if (!DRDU) PetscFunctionReturn(PETSC_SUCCESS);
352814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
353c9ad9de0SHong Zhang   PetscValidHeaderSpecific(U, VEC_CLASSID, 3);
354814e85d6SHong Zhang 
355792fecdfSBarry Smith   PetscCallBack("TS callback DRDU for sensitivity analysis", (*ts->drdufunction)(ts, t, U, DRDU, ts->costintegrandctx));
3563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
357814e85d6SHong Zhang }
358814e85d6SHong Zhang 
359b5b4017aSHong Zhang /*@C
360bcf0153eSBarry Smith   TSComputeDRDPFunction - Deprecated, use `TSGetQuadratureTS()` then `TSComputeRHSJacobianP()`
361814e85d6SHong Zhang 
362d76be551SHong Zhang   Level: deprecated
363814e85d6SHong Zhang 
364814e85d6SHong Zhang @*/
365d71ae5a4SJacob Faibussowitsch PetscErrorCode TSComputeDRDPFunction(TS ts, PetscReal t, Vec U, Vec *DRDP)
366d71ae5a4SJacob Faibussowitsch {
367814e85d6SHong Zhang   PetscFunctionBegin;
3683ba16761SJacob Faibussowitsch   if (!DRDP) PetscFunctionReturn(PETSC_SUCCESS);
369814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
370c9ad9de0SHong Zhang   PetscValidHeaderSpecific(U, VEC_CLASSID, 3);
371814e85d6SHong Zhang 
372792fecdfSBarry Smith   PetscCallBack("TS callback DRDP for sensitivity analysis", (*ts->drdpfunction)(ts, t, U, DRDP, ts->costintegrandctx));
3733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
374814e85d6SHong Zhang }
375814e85d6SHong Zhang 
376b5b4017aSHong Zhang /*@C
377b5b4017aSHong 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.
378b5b4017aSHong Zhang 
379c3339decSBarry Smith   Logically Collective
380b5b4017aSHong Zhang 
381b5b4017aSHong Zhang   Input Parameters:
382bcf0153eSBarry Smith + ts - `TS` context obtained from `TSCreate()`
383b5b4017aSHong Zhang . ihp1 - an array of vectors storing the result of vector-Hessian-vector product for F_UU
384b5b4017aSHong Zhang . hessianproductfunc1 - vector-Hessian-vector product function for F_UU
385b5b4017aSHong Zhang . ihp2 - an array of vectors storing the result of vector-Hessian-vector product for F_UP
386b5b4017aSHong Zhang . hessianproductfunc2 - vector-Hessian-vector product function for F_UP
387b5b4017aSHong Zhang . ihp3 - an array of vectors storing the result of vector-Hessian-vector product for F_PU
388b5b4017aSHong Zhang . hessianproductfunc3 - vector-Hessian-vector product function for F_PU
389b5b4017aSHong Zhang . ihp4 - an array of vectors storing the result of vector-Hessian-vector product for F_PP
390f0fc11ceSJed Brown - hessianproductfunc4 - vector-Hessian-vector product function for F_PP
391b5b4017aSHong Zhang 
392b5b4017aSHong Zhang   Calling sequence of ihessianproductfunc:
393369cf35fSHong Zhang $ ihessianproductfunc (TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx);
394b5b4017aSHong Zhang +   t - current timestep
395b5b4017aSHong Zhang .   U - input vector (current ODE solution)
396369cf35fSHong Zhang .   Vl - an array of input vectors to be left-multiplied with the Hessian
397b5b4017aSHong Zhang .   Vr - input vector to be right-multiplied with the Hessian
398369cf35fSHong Zhang .   VHV - an array of output vectors for vector-Hessian-vector product
399b5b4017aSHong Zhang -   ctx - [optional] user-defined function context
400b5b4017aSHong Zhang 
401b5b4017aSHong Zhang   Level: intermediate
402b5b4017aSHong Zhang 
403369cf35fSHong Zhang   Notes:
404369cf35fSHong Zhang   The first Hessian function and the working array are required.
405369cf35fSHong Zhang   As an example to implement the callback functions, the second callback function calculates the vector-Hessian-vector product
406369cf35fSHong Zhang   $ Vl_n^T*F_UP*Vr
407369cf35fSHong 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.
408369cf35fSHong Zhang   Each entry of F_UP corresponds to the derivative
409369cf35fSHong Zhang   $ F_UP[i][j][k] = \frac{\partial^2 F[i]}{\partial U[j] \partial P[k]}.
410369cf35fSHong 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
411369cf35fSHong Zhang   $ VHV_n[j] = \sum_i \sum_k {Vl_n[i] * F_UP[i][j][k] * Vr[k]}
412369cf35fSHong Zhang   If the cost function is a scalar, there will be only one vector in Vl and VHV.
413b5b4017aSHong Zhang 
414bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`
415b5b4017aSHong Zhang @*/
416d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSetIHessianProduct(TS ts, Vec *ihp1, PetscErrorCode (*ihessianproductfunc1)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *), Vec *ihp2, PetscErrorCode (*ihessianproductfunc2)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *), Vec *ihp3, PetscErrorCode (*ihessianproductfunc3)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *), Vec *ihp4, PetscErrorCode (*ihessianproductfunc4)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *), void *ctx)
417d71ae5a4SJacob Faibussowitsch {
418b5b4017aSHong Zhang   PetscFunctionBegin;
419b5b4017aSHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
420b5b4017aSHong Zhang   PetscValidPointer(ihp1, 2);
421b5b4017aSHong Zhang 
422b5b4017aSHong Zhang   ts->ihessianproductctx = ctx;
423b5b4017aSHong Zhang   if (ihp1) ts->vecs_fuu = ihp1;
424b5b4017aSHong Zhang   if (ihp2) ts->vecs_fup = ihp2;
425b5b4017aSHong Zhang   if (ihp3) ts->vecs_fpu = ihp3;
426b5b4017aSHong Zhang   if (ihp4) ts->vecs_fpp = ihp4;
427b5b4017aSHong Zhang   ts->ihessianproduct_fuu = ihessianproductfunc1;
428b5b4017aSHong Zhang   ts->ihessianproduct_fup = ihessianproductfunc2;
429b5b4017aSHong Zhang   ts->ihessianproduct_fpu = ihessianproductfunc3;
430b5b4017aSHong Zhang   ts->ihessianproduct_fpp = ihessianproductfunc4;
4313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
432b5b4017aSHong Zhang }
433b5b4017aSHong Zhang 
434b5b4017aSHong Zhang /*@C
435b1e111ebSHong Zhang   TSComputeIHessianProductFunctionUU - Runs the user-defined vector-Hessian-vector product function for Fuu.
436b5b4017aSHong Zhang 
437c3339decSBarry Smith   Collective
438b5b4017aSHong Zhang 
439b5b4017aSHong Zhang   Input Parameters:
440bcf0153eSBarry Smith . ts   - The `TS` context obtained from `TSCreate()`
441b5b4017aSHong Zhang 
442b5b4017aSHong Zhang   Level: developer
443b5b4017aSHong Zhang 
444bcf0153eSBarry Smith   Note:
445bcf0153eSBarry Smith   `TSComputeIHessianProductFunctionUU()` is typically used for sensitivity implementation,
446bcf0153eSBarry Smith   so most users would not generally call this routine themselves.
447bcf0153eSBarry Smith 
448bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSSetIHessianProduct()`
449b5b4017aSHong Zhang @*/
450d71ae5a4SJacob Faibussowitsch PetscErrorCode TSComputeIHessianProductFunctionUU(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV)
451d71ae5a4SJacob Faibussowitsch {
452b5b4017aSHong Zhang   PetscFunctionBegin;
4533ba16761SJacob Faibussowitsch   if (!VHV) PetscFunctionReturn(PETSC_SUCCESS);
454b5b4017aSHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
455b5b4017aSHong Zhang   PetscValidHeaderSpecific(U, VEC_CLASSID, 3);
456b5b4017aSHong Zhang 
457792fecdfSBarry Smith   if (ts->ihessianproduct_fuu) PetscCallBack("TS callback IHessianProduct 1 for sensitivity analysis", (*ts->ihessianproduct_fuu)(ts, t, U, Vl, Vr, VHV, ts->ihessianproductctx));
458ef1023bdSBarry Smith 
45967633408SHong Zhang   /* does not consider IMEX for now, so either IHessian or RHSHessian will be calculated, using the same output VHV */
46067633408SHong Zhang   if (ts->rhshessianproduct_guu) {
46167633408SHong Zhang     PetscInt nadj;
4629566063dSJacob Faibussowitsch     PetscCall(TSComputeRHSHessianProductFunctionUU(ts, t, U, Vl, Vr, VHV));
46348a46eb9SPierre Jolivet     for (nadj = 0; nadj < ts->numcost; nadj++) PetscCall(VecScale(VHV[nadj], -1));
46467633408SHong Zhang   }
4653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
466b5b4017aSHong Zhang }
467b5b4017aSHong Zhang 
468b5b4017aSHong Zhang /*@C
469b1e111ebSHong Zhang   TSComputeIHessianProductFunctionUP - Runs the user-defined vector-Hessian-vector product function for Fup.
470b5b4017aSHong Zhang 
471c3339decSBarry Smith   Collective
472b5b4017aSHong Zhang 
473b5b4017aSHong Zhang   Input Parameters:
474bcf0153eSBarry Smith . ts   - The `TS` context obtained from `TSCreate()`
475b5b4017aSHong Zhang 
476b5b4017aSHong Zhang   Level: developer
477b5b4017aSHong Zhang 
478bcf0153eSBarry Smith   Note:
479bcf0153eSBarry Smith   `TSComputeIHessianProductFunctionUP()` is typically used for sensitivity implementation,
480bcf0153eSBarry Smith   so most users would not generally call this routine themselves.
481bcf0153eSBarry Smith 
482bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSSetIHessianProduct()`
483b5b4017aSHong Zhang @*/
484d71ae5a4SJacob Faibussowitsch PetscErrorCode TSComputeIHessianProductFunctionUP(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV)
485d71ae5a4SJacob Faibussowitsch {
486b5b4017aSHong Zhang   PetscFunctionBegin;
4873ba16761SJacob Faibussowitsch   if (!VHV) PetscFunctionReturn(PETSC_SUCCESS);
488b5b4017aSHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
489b5b4017aSHong Zhang   PetscValidHeaderSpecific(U, VEC_CLASSID, 3);
490b5b4017aSHong Zhang 
491792fecdfSBarry Smith   if (ts->ihessianproduct_fup) PetscCallBack("TS callback IHessianProduct 2 for sensitivity analysis", (*ts->ihessianproduct_fup)(ts, t, U, Vl, Vr, VHV, ts->ihessianproductctx));
492ef1023bdSBarry Smith 
49367633408SHong Zhang   /* does not consider IMEX for now, so either IHessian or RHSHessian will be calculated, using the same output VHV */
49467633408SHong Zhang   if (ts->rhshessianproduct_gup) {
49567633408SHong Zhang     PetscInt nadj;
4969566063dSJacob Faibussowitsch     PetscCall(TSComputeRHSHessianProductFunctionUP(ts, t, U, Vl, Vr, VHV));
49748a46eb9SPierre Jolivet     for (nadj = 0; nadj < ts->numcost; nadj++) PetscCall(VecScale(VHV[nadj], -1));
49867633408SHong Zhang   }
4993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
500b5b4017aSHong Zhang }
501b5b4017aSHong Zhang 
502b5b4017aSHong Zhang /*@C
503b1e111ebSHong Zhang   TSComputeIHessianProductFunctionPU - Runs the user-defined vector-Hessian-vector product function for Fpu.
504b5b4017aSHong Zhang 
505c3339decSBarry Smith   Collective
506b5b4017aSHong Zhang 
507b5b4017aSHong Zhang   Input Parameters:
508bcf0153eSBarry Smith . ts   - The `TS` context obtained from `TSCreate()`
509b5b4017aSHong Zhang 
510b5b4017aSHong Zhang   Level: developer
511b5b4017aSHong Zhang 
512bcf0153eSBarry Smith   Note:
513bcf0153eSBarry Smith   `TSComputeIHessianProductFunctionPU()` is typically used for sensitivity implementation,
514bcf0153eSBarry Smith   so most users would not generally call this routine themselves.
515bcf0153eSBarry Smith 
516bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSSetIHessianProduct()`
517b5b4017aSHong Zhang @*/
518d71ae5a4SJacob Faibussowitsch PetscErrorCode TSComputeIHessianProductFunctionPU(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV)
519d71ae5a4SJacob Faibussowitsch {
520b5b4017aSHong Zhang   PetscFunctionBegin;
5213ba16761SJacob Faibussowitsch   if (!VHV) PetscFunctionReturn(PETSC_SUCCESS);
522b5b4017aSHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
523b5b4017aSHong Zhang   PetscValidHeaderSpecific(U, VEC_CLASSID, 3);
524b5b4017aSHong Zhang 
525792fecdfSBarry Smith   if (ts->ihessianproduct_fpu) PetscCallBack("TS callback IHessianProduct 3 for sensitivity analysis", (*ts->ihessianproduct_fpu)(ts, t, U, Vl, Vr, VHV, ts->ihessianproductctx));
526ef1023bdSBarry Smith 
52767633408SHong Zhang   /* does not consider IMEX for now, so either IHessian or RHSHessian will be calculated, using the same output VHV */
52867633408SHong Zhang   if (ts->rhshessianproduct_gpu) {
52967633408SHong Zhang     PetscInt nadj;
5309566063dSJacob Faibussowitsch     PetscCall(TSComputeRHSHessianProductFunctionPU(ts, t, U, Vl, Vr, VHV));
53148a46eb9SPierre Jolivet     for (nadj = 0; nadj < ts->numcost; nadj++) PetscCall(VecScale(VHV[nadj], -1));
53267633408SHong Zhang   }
5333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
534b5b4017aSHong Zhang }
535b5b4017aSHong Zhang 
536b5b4017aSHong Zhang /*@C
537b1e111ebSHong Zhang   TSComputeIHessianProductFunctionPP - Runs the user-defined vector-Hessian-vector product function for Fpp.
538b5b4017aSHong Zhang 
539c3339decSBarry Smith   Collective
540b5b4017aSHong Zhang 
541b5b4017aSHong Zhang   Input Parameters:
542bcf0153eSBarry Smith . ts   - The `TS` context obtained from `TSCreate()`
543b5b4017aSHong Zhang 
544b5b4017aSHong Zhang   Level: developer
545b5b4017aSHong Zhang 
546bcf0153eSBarry Smith   Note:
547bcf0153eSBarry Smith   `TSComputeIHessianProductFunctionPP()` is typically used for sensitivity implementation,
548bcf0153eSBarry Smith   so most users would not generally call this routine themselves.
549bcf0153eSBarry Smith 
550bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSSetIHessianProduct()`
551b5b4017aSHong Zhang @*/
552d71ae5a4SJacob Faibussowitsch PetscErrorCode TSComputeIHessianProductFunctionPP(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV)
553d71ae5a4SJacob Faibussowitsch {
554b5b4017aSHong Zhang   PetscFunctionBegin;
5553ba16761SJacob Faibussowitsch   if (!VHV) PetscFunctionReturn(PETSC_SUCCESS);
556b5b4017aSHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
557b5b4017aSHong Zhang   PetscValidHeaderSpecific(U, VEC_CLASSID, 3);
558b5b4017aSHong Zhang 
559792fecdfSBarry Smith   if (ts->ihessianproduct_fpp) PetscCallBack("TS callback IHessianProduct 3 for sensitivity analysis", (*ts->ihessianproduct_fpp)(ts, t, U, Vl, Vr, VHV, ts->ihessianproductctx));
560ef1023bdSBarry Smith 
56167633408SHong Zhang   /* does not consider IMEX for now, so either IHessian or RHSHessian will be calculated, using the same output VHV */
56267633408SHong Zhang   if (ts->rhshessianproduct_gpp) {
56367633408SHong Zhang     PetscInt nadj;
5649566063dSJacob Faibussowitsch     PetscCall(TSComputeRHSHessianProductFunctionPP(ts, t, U, Vl, Vr, VHV));
56548a46eb9SPierre Jolivet     for (nadj = 0; nadj < ts->numcost; nadj++) PetscCall(VecScale(VHV[nadj], -1));
56667633408SHong Zhang   }
5673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
568b5b4017aSHong Zhang }
569b5b4017aSHong Zhang 
57013af1a74SHong Zhang /*@C
5714c500f23SPierre 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.
57213af1a74SHong Zhang 
573c3339decSBarry Smith   Logically Collective
57413af1a74SHong Zhang 
57513af1a74SHong Zhang   Input Parameters:
576bcf0153eSBarry Smith + ts - `TS` context obtained from `TSCreate()`
57713af1a74SHong Zhang . rhshp1 - an array of vectors storing the result of vector-Hessian-vector product for G_UU
57813af1a74SHong Zhang . hessianproductfunc1 - vector-Hessian-vector product function for G_UU
57913af1a74SHong Zhang . rhshp2 - an array of vectors storing the result of vector-Hessian-vector product for G_UP
58013af1a74SHong Zhang . hessianproductfunc2 - vector-Hessian-vector product function for G_UP
58113af1a74SHong Zhang . rhshp3 - an array of vectors storing the result of vector-Hessian-vector product for G_PU
58213af1a74SHong Zhang . hessianproductfunc3 - vector-Hessian-vector product function for G_PU
58313af1a74SHong Zhang . rhshp4 - an array of vectors storing the result of vector-Hessian-vector product for G_PP
584f0fc11ceSJed Brown - hessianproductfunc4 - vector-Hessian-vector product function for G_PP
58513af1a74SHong Zhang 
58613af1a74SHong Zhang   Calling sequence of ihessianproductfunc:
587369cf35fSHong Zhang $ rhshessianproductfunc (TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx);
58813af1a74SHong Zhang +   t - current timestep
58913af1a74SHong Zhang .   U - input vector (current ODE solution)
590369cf35fSHong Zhang .   Vl - an array of input vectors to be left-multiplied with the Hessian
59113af1a74SHong Zhang .   Vr - input vector to be right-multiplied with the Hessian
592369cf35fSHong Zhang .   VHV - an array of output vectors for vector-Hessian-vector product
59313af1a74SHong Zhang -   ctx - [optional] user-defined function context
59413af1a74SHong Zhang 
59513af1a74SHong Zhang   Level: intermediate
59613af1a74SHong Zhang 
597369cf35fSHong Zhang   Notes:
598369cf35fSHong Zhang   The first Hessian function and the working array are required.
599369cf35fSHong Zhang   As an example to implement the callback functions, the second callback function calculates the vector-Hessian-vector product
600369cf35fSHong Zhang   $ Vl_n^T*G_UP*Vr
601369cf35fSHong 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.
602369cf35fSHong Zhang   Each entry of G_UP corresponds to the derivative
603369cf35fSHong Zhang   $ G_UP[i][j][k] = \frac{\partial^2 G[i]}{\partial U[j] \partial P[k]}.
604369cf35fSHong 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
605369cf35fSHong Zhang   $ VHV_n[j] = \sum_i \sum_k {Vl_n[i] * G_UP[i][j][k] * Vr[k]}
606369cf35fSHong Zhang   If the cost function is a scalar, there will be only one vector in Vl and VHV.
60713af1a74SHong Zhang 
608bcf0153eSBarry Smith .seealso: `TS`
60913af1a74SHong Zhang @*/
610d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSetRHSHessianProduct(TS ts, Vec *rhshp1, PetscErrorCode (*rhshessianproductfunc1)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *), Vec *rhshp2, PetscErrorCode (*rhshessianproductfunc2)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *), Vec *rhshp3, PetscErrorCode (*rhshessianproductfunc3)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *), Vec *rhshp4, PetscErrorCode (*rhshessianproductfunc4)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *), void *ctx)
611d71ae5a4SJacob Faibussowitsch {
61213af1a74SHong Zhang   PetscFunctionBegin;
61313af1a74SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
61413af1a74SHong Zhang   PetscValidPointer(rhshp1, 2);
61513af1a74SHong Zhang 
61613af1a74SHong Zhang   ts->rhshessianproductctx = ctx;
61713af1a74SHong Zhang   if (rhshp1) ts->vecs_guu = rhshp1;
61813af1a74SHong Zhang   if (rhshp2) ts->vecs_gup = rhshp2;
61913af1a74SHong Zhang   if (rhshp3) ts->vecs_gpu = rhshp3;
62013af1a74SHong Zhang   if (rhshp4) ts->vecs_gpp = rhshp4;
62113af1a74SHong Zhang   ts->rhshessianproduct_guu = rhshessianproductfunc1;
62213af1a74SHong Zhang   ts->rhshessianproduct_gup = rhshessianproductfunc2;
62313af1a74SHong Zhang   ts->rhshessianproduct_gpu = rhshessianproductfunc3;
62413af1a74SHong Zhang   ts->rhshessianproduct_gpp = rhshessianproductfunc4;
6253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
62613af1a74SHong Zhang }
62713af1a74SHong Zhang 
62813af1a74SHong Zhang /*@C
629b1e111ebSHong Zhang   TSComputeRHSHessianProductFunctionUU - Runs the user-defined vector-Hessian-vector product function for Guu.
63013af1a74SHong Zhang 
631c3339decSBarry Smith   Collective
63213af1a74SHong Zhang 
63313af1a74SHong Zhang   Input Parameters:
634bcf0153eSBarry Smith . ts   - The `TS` context obtained from `TSCreate()`
63513af1a74SHong Zhang 
63613af1a74SHong Zhang   Level: developer
63713af1a74SHong Zhang 
638bcf0153eSBarry Smith   Note:
639bcf0153eSBarry Smith   `TSComputeRHSHessianProductFunctionUU()` is typically used for sensitivity implementation,
640bcf0153eSBarry Smith   so most users would not generally call this routine themselves.
641bcf0153eSBarry Smith 
642bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSSetRHSHessianProduct()`
64313af1a74SHong Zhang @*/
644d71ae5a4SJacob Faibussowitsch PetscErrorCode TSComputeRHSHessianProductFunctionUU(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV)
645d71ae5a4SJacob Faibussowitsch {
64613af1a74SHong Zhang   PetscFunctionBegin;
6473ba16761SJacob Faibussowitsch   if (!VHV) PetscFunctionReturn(PETSC_SUCCESS);
64813af1a74SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
64913af1a74SHong Zhang   PetscValidHeaderSpecific(U, VEC_CLASSID, 3);
65013af1a74SHong Zhang 
651792fecdfSBarry Smith   PetscCallBack("TS callback RHSHessianProduct 1 for sensitivity analysis", (*ts->rhshessianproduct_guu)(ts, t, U, Vl, Vr, VHV, ts->rhshessianproductctx));
6523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
65313af1a74SHong Zhang }
65413af1a74SHong Zhang 
65513af1a74SHong Zhang /*@C
656b1e111ebSHong Zhang   TSComputeRHSHessianProductFunctionUP - Runs the user-defined vector-Hessian-vector product function for Gup.
65713af1a74SHong Zhang 
658c3339decSBarry Smith   Collective
65913af1a74SHong Zhang 
66013af1a74SHong Zhang   Input Parameters:
661bcf0153eSBarry Smith . ts   - The `TS` context obtained from `TSCreate()`
66213af1a74SHong Zhang 
66313af1a74SHong Zhang   Level: developer
66413af1a74SHong Zhang 
665bcf0153eSBarry Smith   Note:
666bcf0153eSBarry Smith   `TSComputeRHSHessianProductFunctionUP()` is typically used for sensitivity implementation,
667bcf0153eSBarry Smith   so most users would not generally call this routine themselves.
668bcf0153eSBarry Smith 
669bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSSetRHSHessianProduct()`
67013af1a74SHong Zhang @*/
671d71ae5a4SJacob Faibussowitsch PetscErrorCode TSComputeRHSHessianProductFunctionUP(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV)
672d71ae5a4SJacob Faibussowitsch {
67313af1a74SHong Zhang   PetscFunctionBegin;
6743ba16761SJacob Faibussowitsch   if (!VHV) PetscFunctionReturn(PETSC_SUCCESS);
67513af1a74SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
67613af1a74SHong Zhang   PetscValidHeaderSpecific(U, VEC_CLASSID, 3);
67713af1a74SHong Zhang 
678792fecdfSBarry Smith   PetscCallBack("TS callback RHSHessianProduct 2 for sensitivity analysis", (*ts->rhshessianproduct_gup)(ts, t, U, Vl, Vr, VHV, ts->rhshessianproductctx));
6793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
68013af1a74SHong Zhang }
68113af1a74SHong Zhang 
68213af1a74SHong Zhang /*@C
683b1e111ebSHong Zhang   TSComputeRHSHessianProductFunctionPU - Runs the user-defined vector-Hessian-vector product function for Gpu.
68413af1a74SHong Zhang 
685c3339decSBarry Smith   Collective
68613af1a74SHong Zhang 
68713af1a74SHong Zhang   Input Parameters:
688bcf0153eSBarry Smith . ts   - The `TS` context obtained from `TSCreate()`
68913af1a74SHong Zhang 
69013af1a74SHong Zhang   Level: developer
69113af1a74SHong Zhang 
692bcf0153eSBarry Smith   Note:
693bcf0153eSBarry Smith   `TSComputeRHSHessianProductFunctionPU()` is typically used for sensitivity implementation,
694bcf0153eSBarry Smith   so most users would not generally call this routine themselves.
695bcf0153eSBarry Smith 
696bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSSetRHSHessianProduct()`
69713af1a74SHong Zhang @*/
698d71ae5a4SJacob Faibussowitsch PetscErrorCode TSComputeRHSHessianProductFunctionPU(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV)
699d71ae5a4SJacob Faibussowitsch {
70013af1a74SHong Zhang   PetscFunctionBegin;
7013ba16761SJacob Faibussowitsch   if (!VHV) PetscFunctionReturn(PETSC_SUCCESS);
70213af1a74SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
70313af1a74SHong Zhang   PetscValidHeaderSpecific(U, VEC_CLASSID, 3);
70413af1a74SHong Zhang 
705792fecdfSBarry Smith   PetscCallBack("TS callback RHSHessianProduct 3 for sensitivity analysis", (*ts->rhshessianproduct_gpu)(ts, t, U, Vl, Vr, VHV, ts->rhshessianproductctx));
7063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
70713af1a74SHong Zhang }
70813af1a74SHong Zhang 
70913af1a74SHong Zhang /*@C
710b1e111ebSHong Zhang   TSComputeRHSHessianProductFunctionPP - Runs the user-defined vector-Hessian-vector product function for Gpp.
71113af1a74SHong Zhang 
712c3339decSBarry Smith   Collective
71313af1a74SHong Zhang 
71413af1a74SHong Zhang   Input Parameters:
715bcf0153eSBarry Smith . ts   - The `TS` context obtained from `TSCreate()`
71613af1a74SHong Zhang 
71713af1a74SHong Zhang   Level: developer
71813af1a74SHong Zhang 
719bcf0153eSBarry Smith   Note:
720bcf0153eSBarry Smith   `TSComputeRHSHessianProductFunctionPP()` is typically used for sensitivity implementation,
721bcf0153eSBarry Smith   so most users would not generally call this routine themselves.
722bcf0153eSBarry Smith 
723bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSSetRHSHessianProduct()`
72413af1a74SHong Zhang @*/
725d71ae5a4SJacob Faibussowitsch PetscErrorCode TSComputeRHSHessianProductFunctionPP(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV)
726d71ae5a4SJacob Faibussowitsch {
72713af1a74SHong Zhang   PetscFunctionBegin;
7283ba16761SJacob Faibussowitsch   if (!VHV) PetscFunctionReturn(PETSC_SUCCESS);
72913af1a74SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
73013af1a74SHong Zhang   PetscValidHeaderSpecific(U, VEC_CLASSID, 3);
73113af1a74SHong Zhang 
732792fecdfSBarry Smith   PetscCallBack("TS callback RHSHessianProduct 3 for sensitivity analysis", (*ts->rhshessianproduct_gpp)(ts, t, U, Vl, Vr, VHV, ts->rhshessianproductctx));
7333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
73413af1a74SHong Zhang }
73513af1a74SHong Zhang 
736814e85d6SHong Zhang /* --------------------------- Adjoint sensitivity ---------------------------*/
737814e85d6SHong Zhang 
738814e85d6SHong Zhang /*@
739814e85d6SHong 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
740bcf0153eSBarry Smith       for use by the `TS` adjoint routines.
741814e85d6SHong Zhang 
742c3339decSBarry Smith    Logically Collective
743814e85d6SHong Zhang 
744814e85d6SHong Zhang    Input Parameters:
745bcf0153eSBarry Smith +  ts - the `TS` context obtained from `TSCreate()`
746d2fd7bfcSHong Zhang .  numcost - number of gradients to be computed, this is the number of cost functions
747814e85d6SHong 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
748814e85d6SHong Zhang -  mu - gradients with respect to the parameters, the number of entries in these vectors is the same as the number of parameters
749814e85d6SHong Zhang 
750814e85d6SHong Zhang    Level: beginner
751814e85d6SHong Zhang 
752814e85d6SHong Zhang    Notes:
75335cb6cd3SPierre Jolivet     the entries in these vectors must be correctly initialized with the values lambda_i = df/dy|finaltime  mu_i = df/dp|finaltime
754814e85d6SHong Zhang 
75535cb6cd3SPierre Jolivet    After `TSAdjointSolve()` is called the lambda and the mu contain the computed sensitivities
756814e85d6SHong Zhang 
757bcf0153eSBarry Smith .seealso: `TS`, `TSAdjointSolve()`, `TSGetCostGradients()`
758814e85d6SHong Zhang @*/
759d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSetCostGradients(TS ts, PetscInt numcost, Vec *lambda, Vec *mu)
760d71ae5a4SJacob Faibussowitsch {
761814e85d6SHong Zhang   PetscFunctionBegin;
762814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
763064a246eSJacob Faibussowitsch   PetscValidPointer(lambda, 3);
764814e85d6SHong Zhang   ts->vecs_sensi  = lambda;
765814e85d6SHong Zhang   ts->vecs_sensip = mu;
7663c633725SBarry Smith   PetscCheck(!ts->numcost || ts->numcost == numcost, PetscObjectComm((PetscObject)ts), PETSC_ERR_USER, "The number of cost functions (2nd parameter of TSSetCostIntegrand()) is inconsistent with the one set by TSSetCostIntegrand");
767814e85d6SHong Zhang   ts->numcost = numcost;
7683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
769814e85d6SHong Zhang }
770814e85d6SHong Zhang 
771814e85d6SHong Zhang /*@
772bcf0153eSBarry Smith    TSGetCostGradients - Returns the gradients from the `TSAdjointSolve()`
773814e85d6SHong Zhang 
774bcf0153eSBarry Smith    Not Collective, but the vectors returned are parallel if `TS` is parallel
775814e85d6SHong Zhang 
776814e85d6SHong Zhang    Input Parameter:
777bcf0153eSBarry Smith .  ts - the `TS` context obtained from `TSCreate()`
778814e85d6SHong Zhang 
779d8d19677SJose E. Roman    Output Parameters:
7806b867d5aSJose E. Roman +  numcost - size of returned arrays
7816b867d5aSJose E. Roman .  lambda - vectors containing the gradients of the cost functions with respect to the ODE/DAE solution variables
782814e85d6SHong Zhang -  mu - vectors containing the gradients of the cost functions with respect to the problem parameters
783814e85d6SHong Zhang 
784814e85d6SHong Zhang    Level: intermediate
785814e85d6SHong Zhang 
786bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSAdjointSolve()`, `TSSetCostGradients()`
787814e85d6SHong Zhang @*/
788d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGetCostGradients(TS ts, PetscInt *numcost, Vec **lambda, Vec **mu)
789d71ae5a4SJacob Faibussowitsch {
790814e85d6SHong Zhang   PetscFunctionBegin;
791814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
792814e85d6SHong Zhang   if (numcost) *numcost = ts->numcost;
793814e85d6SHong Zhang   if (lambda) *lambda = ts->vecs_sensi;
794814e85d6SHong Zhang   if (mu) *mu = ts->vecs_sensip;
7953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
796814e85d6SHong Zhang }
797814e85d6SHong Zhang 
798814e85d6SHong Zhang /*@
799b5b4017aSHong 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
800bcf0153eSBarry Smith    for use by the `TS` adjoint routines.
801b5b4017aSHong Zhang 
802c3339decSBarry Smith    Logically Collective
803b5b4017aSHong Zhang 
804b5b4017aSHong Zhang    Input Parameters:
805bcf0153eSBarry Smith +  ts - the `TS` context obtained from `TSCreate()`
806b5b4017aSHong Zhang .  numcost - number of cost functions
807b5b4017aSHong 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
808b5b4017aSHong 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
809b5b4017aSHong Zhang -  dir - the direction vector that are multiplied with the Hessian of the cost functions
810b5b4017aSHong Zhang 
811b5b4017aSHong Zhang    Level: beginner
812b5b4017aSHong Zhang 
813bcf0153eSBarry Smith    Notes:
814bcf0153eSBarry Smith    Hessian of the cost function is completely different from Hessian of the ODE/DAE system
815b5b4017aSHong Zhang 
816bcf0153eSBarry Smith    For second-order adjoint, one needs to call this function and then `TSAdjointSetForward()` before `TSSolve()`.
817b5b4017aSHong Zhang 
81835cb6cd3SPierre Jolivet    After `TSAdjointSolve()` is called, the lambda2 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.
819b5b4017aSHong Zhang 
8203fca9621SHong Zhang    Passing NULL for lambda2 disables the second-order calculation.
821bcf0153eSBarry Smith 
822bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSAdjointSolve()`, `TSAdjointSetForward()`
823b5b4017aSHong Zhang @*/
824d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSetCostHessianProducts(TS ts, PetscInt numcost, Vec *lambda2, Vec *mu2, Vec dir)
825d71ae5a4SJacob Faibussowitsch {
826b5b4017aSHong Zhang   PetscFunctionBegin;
827b5b4017aSHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
8283c633725SBarry Smith   PetscCheck(!ts->numcost || ts->numcost == numcost, PetscObjectComm((PetscObject)ts), PETSC_ERR_USER, "The number of cost functions (2nd parameter of TSSetCostIntegrand()) is inconsistent with the one set by TSSetCostIntegrand");
829b5b4017aSHong Zhang   ts->numcost      = numcost;
830b5b4017aSHong Zhang   ts->vecs_sensi2  = lambda2;
83133f72c5dSHong Zhang   ts->vecs_sensi2p = mu2;
832b5b4017aSHong Zhang   ts->vec_dir      = dir;
8333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
834b5b4017aSHong Zhang }
835b5b4017aSHong Zhang 
836b5b4017aSHong Zhang /*@
837bcf0153eSBarry Smith    TSGetCostHessianProducts - Returns the gradients from the `TSAdjointSolve()`
838b5b4017aSHong Zhang 
839bcf0153eSBarry Smith    Not Collective, but vectors returned are parallel if `TS` is parallel
840b5b4017aSHong Zhang 
841b5b4017aSHong Zhang    Input Parameter:
842bcf0153eSBarry Smith .  ts - the `TS` context obtained from `TSCreate()`
843b5b4017aSHong Zhang 
844d8d19677SJose E. Roman    Output Parameters:
845b5b4017aSHong Zhang +  numcost - number of cost functions
846b5b4017aSHong 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
847b5b4017aSHong 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
848b5b4017aSHong Zhang -  dir - the direction vector that are multiplied with the Hessian of the cost functions
849b5b4017aSHong Zhang 
850b5b4017aSHong Zhang    Level: intermediate
851b5b4017aSHong Zhang 
852bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSAdjointSolve()`, `TSSetCostHessianProducts()`
853b5b4017aSHong Zhang @*/
854d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGetCostHessianProducts(TS ts, PetscInt *numcost, Vec **lambda2, Vec **mu2, Vec *dir)
855d71ae5a4SJacob Faibussowitsch {
856b5b4017aSHong Zhang   PetscFunctionBegin;
857b5b4017aSHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
858b5b4017aSHong Zhang   if (numcost) *numcost = ts->numcost;
859b5b4017aSHong Zhang   if (lambda2) *lambda2 = ts->vecs_sensi2;
86033f72c5dSHong Zhang   if (mu2) *mu2 = ts->vecs_sensi2p;
861b5b4017aSHong Zhang   if (dir) *dir = ts->vec_dir;
8623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
863b5b4017aSHong Zhang }
864b5b4017aSHong Zhang 
865b5b4017aSHong Zhang /*@
866ecf68647SHong Zhang   TSAdjointSetForward - Trigger the tangent linear solver and initialize the forward sensitivities
867b5b4017aSHong Zhang 
868c3339decSBarry Smith   Logically Collective
869b5b4017aSHong Zhang 
870b5b4017aSHong Zhang   Input Parameters:
871bcf0153eSBarry Smith +  ts - the `TS` context obtained from `TSCreate()`
872b5b4017aSHong Zhang -  didp - the derivative of initial values w.r.t. parameters
873b5b4017aSHong Zhang 
874b5b4017aSHong Zhang   Level: intermediate
875b5b4017aSHong Zhang 
876bcf0153eSBarry Smith   Notes:
877bcf0153eSBarry Smith   When computing sensitivies w.r.t. initial condition, set didp to NULL so that the solver will take it as an identity matrix mathematically.
878bcf0153eSBarry Smith   `TS` adjoint does not reset the tangent linear solver automatically, `TSAdjointResetForward()` should be called to reset the tangent linear solver.
879b5b4017aSHong Zhang 
880bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSAdjointSolve()`, `TSSetCostHessianProducts()`, `TSAdjointResetForward()`
881b5b4017aSHong Zhang @*/
882d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointSetForward(TS ts, Mat didp)
883d71ae5a4SJacob Faibussowitsch {
88433f72c5dSHong Zhang   Mat          A;
88533f72c5dSHong Zhang   Vec          sp;
88633f72c5dSHong Zhang   PetscScalar *xarr;
88733f72c5dSHong Zhang   PetscInt     lsize;
888b5b4017aSHong Zhang 
889b5b4017aSHong Zhang   PetscFunctionBegin;
890b5b4017aSHong Zhang   ts->forward_solve = PETSC_TRUE; /* turn on tangent linear mode */
8913c633725SBarry Smith   PetscCheck(ts->vecs_sensi2, PetscObjectComm((PetscObject)ts), PETSC_ERR_USER, "Must call TSSetCostHessianProducts() first");
8923c633725SBarry Smith   PetscCheck(ts->vec_dir, PetscObjectComm((PetscObject)ts), PETSC_ERR_USER, "Directional vector is missing. Call TSSetCostHessianProducts() to set it.");
89333f72c5dSHong Zhang   /* create a single-column dense matrix */
8949566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(ts->vec_sol, &lsize));
8959566063dSJacob Faibussowitsch   PetscCall(MatCreateDense(PetscObjectComm((PetscObject)ts), lsize, PETSC_DECIDE, PETSC_DECIDE, 1, NULL, &A));
89633f72c5dSHong Zhang 
8979566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &sp));
8989566063dSJacob Faibussowitsch   PetscCall(MatDenseGetColumn(A, 0, &xarr));
8999566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(sp, xarr));
900ecf68647SHong Zhang   if (ts->vecs_sensi2p) { /* tangent linear variable initialized as 2*dIdP*dir */
90133f72c5dSHong Zhang     if (didp) {
9029566063dSJacob Faibussowitsch       PetscCall(MatMult(didp, ts->vec_dir, sp));
9039566063dSJacob Faibussowitsch       PetscCall(VecScale(sp, 2.));
90433f72c5dSHong Zhang     } else {
9059566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(sp));
90633f72c5dSHong Zhang     }
907ecf68647SHong Zhang   } else { /* tangent linear variable initialized as dir */
9089566063dSJacob Faibussowitsch     PetscCall(VecCopy(ts->vec_dir, sp));
90933f72c5dSHong Zhang   }
9109566063dSJacob Faibussowitsch   PetscCall(VecResetArray(sp));
9119566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreColumn(A, &xarr));
9129566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&sp));
91333f72c5dSHong Zhang 
9149566063dSJacob Faibussowitsch   PetscCall(TSForwardSetInitialSensitivities(ts, A)); /* if didp is NULL, identity matrix is assumed */
91533f72c5dSHong Zhang 
9169566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
9173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
918b5b4017aSHong Zhang }
919b5b4017aSHong Zhang 
920b5b4017aSHong Zhang /*@
921ecf68647SHong Zhang   TSAdjointResetForward - Reset the tangent linear solver and destroy the tangent linear context
922ecf68647SHong Zhang 
923c3339decSBarry Smith   Logically Collective
924ecf68647SHong Zhang 
925ecf68647SHong Zhang   Input Parameters:
926bcf0153eSBarry Smith .  ts - the `TS` context obtained from `TSCreate()`
927ecf68647SHong Zhang 
928ecf68647SHong Zhang   Level: intermediate
929ecf68647SHong Zhang 
930bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSAdjointSetForward()`
931ecf68647SHong Zhang @*/
932d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointResetForward(TS ts)
933d71ae5a4SJacob Faibussowitsch {
934ecf68647SHong Zhang   PetscFunctionBegin;
935ecf68647SHong Zhang   ts->forward_solve = PETSC_FALSE; /* turn off tangent linear mode */
9369566063dSJacob Faibussowitsch   PetscCall(TSForwardReset(ts));
9373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
938ecf68647SHong Zhang }
939ecf68647SHong Zhang 
940ecf68647SHong Zhang /*@
941814e85d6SHong Zhang    TSAdjointSetUp - Sets up the internal data structures for the later use
942814e85d6SHong Zhang    of an adjoint solver
943814e85d6SHong Zhang 
944c3339decSBarry Smith    Collective
945814e85d6SHong Zhang 
946814e85d6SHong Zhang    Input Parameter:
947bcf0153eSBarry Smith .  ts - the `TS` context obtained from `TSCreate()`
948814e85d6SHong Zhang 
949814e85d6SHong Zhang    Level: advanced
950814e85d6SHong Zhang 
951bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSCreate()`, `TSAdjointStep()`, `TSSetCostGradients()`
952814e85d6SHong Zhang @*/
953d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointSetUp(TS ts)
954d71ae5a4SJacob Faibussowitsch {
955881c1a9bSHong Zhang   TSTrajectory tj;
956881c1a9bSHong Zhang   PetscBool    match;
957814e85d6SHong Zhang 
958814e85d6SHong Zhang   PetscFunctionBegin;
959814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
9603ba16761SJacob Faibussowitsch   if (ts->adjointsetupcalled) PetscFunctionReturn(PETSC_SUCCESS);
9613c633725SBarry Smith   PetscCheck(ts->vecs_sensi, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_WRONGSTATE, "Must call TSSetCostGradients() first");
9623c633725SBarry Smith   PetscCheck(!ts->vecs_sensip || ts->Jacp || ts->Jacprhs, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_WRONGSTATE, "Must call TSSetRHSJacobianP() or TSSetIJacobianP() first");
9639566063dSJacob Faibussowitsch   PetscCall(TSGetTrajectory(ts, &tj));
9649566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)tj, TSTRAJECTORYBASIC, &match));
965881c1a9bSHong Zhang   if (match) {
966881c1a9bSHong Zhang     PetscBool solution_only;
9679566063dSJacob Faibussowitsch     PetscCall(TSTrajectoryGetSolutionOnly(tj, &solution_only));
9683c633725SBarry Smith     PetscCheck(!solution_only, 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");
969881c1a9bSHong Zhang   }
9709566063dSJacob Faibussowitsch   PetscCall(TSTrajectorySetUseHistory(tj, PETSC_FALSE)); /* not use TSHistory */
97133f72c5dSHong Zhang 
972cd4cee2dSHong Zhang   if (ts->quadraturets) { /* if there is integral in the cost function */
9739566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(ts->vecs_sensi[0], &ts->vec_drdu_col));
97448a46eb9SPierre Jolivet     if (ts->vecs_sensip) PetscCall(VecDuplicate(ts->vecs_sensip[0], &ts->vec_drdp_col));
975814e85d6SHong Zhang   }
976814e85d6SHong Zhang 
977dbbe0bcdSBarry Smith   PetscTryTypeMethod(ts, adjointsetup);
978814e85d6SHong Zhang   ts->adjointsetupcalled = PETSC_TRUE;
9793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
980814e85d6SHong Zhang }
981814e85d6SHong Zhang 
982814e85d6SHong Zhang /*@
983bcf0153eSBarry Smith   TSAdjointReset - Resets a `TS` adjoint context and removes any allocated `Vec`s and `Mat`s.
984ece46509SHong Zhang 
985c3339decSBarry Smith    Collective
986ece46509SHong Zhang 
987ece46509SHong Zhang    Input Parameter:
988bcf0153eSBarry Smith .  ts - the `TS` context obtained from `TSCreate()`
989ece46509SHong Zhang 
990ece46509SHong Zhang    Level: beginner
991ece46509SHong Zhang 
992bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSCreate()`, `TSAdjointSetUp()`, `TSADestroy()`
993ece46509SHong Zhang @*/
994d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointReset(TS ts)
995d71ae5a4SJacob Faibussowitsch {
996ece46509SHong Zhang   PetscFunctionBegin;
997ece46509SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
998dbbe0bcdSBarry Smith   PetscTryTypeMethod(ts, adjointreset);
9997207e0fdSHong Zhang   if (ts->quadraturets) { /* if there is integral in the cost function */
10009566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&ts->vec_drdu_col));
100148a46eb9SPierre Jolivet     if (ts->vecs_sensip) PetscCall(VecDestroy(&ts->vec_drdp_col));
10027207e0fdSHong Zhang   }
1003ece46509SHong Zhang   ts->vecs_sensi         = NULL;
1004ece46509SHong Zhang   ts->vecs_sensip        = NULL;
1005ece46509SHong Zhang   ts->vecs_sensi2        = NULL;
100633f72c5dSHong Zhang   ts->vecs_sensi2p       = NULL;
1007ece46509SHong Zhang   ts->vec_dir            = NULL;
1008ece46509SHong Zhang   ts->adjointsetupcalled = PETSC_FALSE;
10093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1010ece46509SHong Zhang }
1011ece46509SHong Zhang 
1012ece46509SHong Zhang /*@
1013814e85d6SHong Zhang    TSAdjointSetSteps - Sets the number of steps the adjoint solver should take backward in time
1014814e85d6SHong Zhang 
1015c3339decSBarry Smith    Logically Collective
1016814e85d6SHong Zhang 
1017814e85d6SHong Zhang    Input Parameters:
1018bcf0153eSBarry Smith +  ts - the `TS` context obtained from `TSCreate()`
1019a2b725a8SWilliam Gropp -  steps - number of steps to use
1020814e85d6SHong Zhang 
1021814e85d6SHong Zhang    Level: intermediate
1022814e85d6SHong Zhang 
1023814e85d6SHong Zhang    Notes:
1024bcf0153eSBarry Smith     Normally one does not call this and `TSAdjointSolve()` integrates back to the original timestep. One can call this
1025814e85d6SHong Zhang           so as to integrate back to less than the original timestep
1026814e85d6SHong Zhang 
1027bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSAdjointSolve()`, `TS`, `TSSetExactFinalTime()`
1028814e85d6SHong Zhang @*/
1029d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointSetSteps(TS ts, PetscInt steps)
1030d71ae5a4SJacob Faibussowitsch {
1031814e85d6SHong Zhang   PetscFunctionBegin;
1032814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1033814e85d6SHong Zhang   PetscValidLogicalCollectiveInt(ts, steps, 2);
10343c633725SBarry Smith   PetscCheck(steps >= 0, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_OUTOFRANGE, "Cannot step back a negative number of steps");
10353c633725SBarry Smith   PetscCheck(steps <= ts->steps, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_OUTOFRANGE, "Cannot step back more than the total number of forward steps");
1036814e85d6SHong Zhang   ts->adjoint_max_steps = steps;
10373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1038814e85d6SHong Zhang }
1039814e85d6SHong Zhang 
1040814e85d6SHong Zhang /*@C
1041bcf0153eSBarry Smith   TSAdjointSetRHSJacobian - Deprecated, use `TSSetRHSJacobianP()`
1042814e85d6SHong Zhang 
1043814e85d6SHong Zhang   Level: deprecated
1044814e85d6SHong Zhang 
1045814e85d6SHong Zhang @*/
1046d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointSetRHSJacobian(TS ts, Mat Amat, PetscErrorCode (*func)(TS, PetscReal, Vec, Mat, void *), void *ctx)
1047d71ae5a4SJacob Faibussowitsch {
1048814e85d6SHong Zhang   PetscFunctionBegin;
1049814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1050814e85d6SHong Zhang   PetscValidHeaderSpecific(Amat, MAT_CLASSID, 2);
1051814e85d6SHong Zhang 
1052814e85d6SHong Zhang   ts->rhsjacobianp    = func;
1053814e85d6SHong Zhang   ts->rhsjacobianpctx = ctx;
1054814e85d6SHong Zhang   if (Amat) {
10559566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)Amat));
10569566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&ts->Jacp));
1057814e85d6SHong Zhang     ts->Jacp = Amat;
1058814e85d6SHong Zhang   }
10593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1060814e85d6SHong Zhang }
1061814e85d6SHong Zhang 
1062814e85d6SHong Zhang /*@C
1063bcf0153eSBarry Smith   TSAdjointComputeRHSJacobian - Deprecated, use `TSComputeRHSJacobianP()`
1064814e85d6SHong Zhang 
1065814e85d6SHong Zhang   Level: deprecated
1066814e85d6SHong Zhang 
1067814e85d6SHong Zhang @*/
1068d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointComputeRHSJacobian(TS ts, PetscReal t, Vec U, Mat Amat)
1069d71ae5a4SJacob Faibussowitsch {
1070814e85d6SHong Zhang   PetscFunctionBegin;
1071814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1072c9ad9de0SHong Zhang   PetscValidHeaderSpecific(U, VEC_CLASSID, 3);
1073814e85d6SHong Zhang   PetscValidPointer(Amat, 4);
1074814e85d6SHong Zhang 
1075792fecdfSBarry Smith   PetscCallBack("TS callback JacobianP for sensitivity analysis", (*ts->rhsjacobianp)(ts, t, U, Amat, ts->rhsjacobianpctx));
10763ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1077814e85d6SHong Zhang }
1078814e85d6SHong Zhang 
1079814e85d6SHong Zhang /*@
1080bcf0153eSBarry Smith   TSAdjointComputeDRDYFunction - Deprecated, use `TSGetQuadratureTS()` then `TSComputeRHSJacobian()`
1081814e85d6SHong Zhang 
1082814e85d6SHong Zhang   Level: deprecated
1083814e85d6SHong Zhang 
1084814e85d6SHong Zhang @*/
1085d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointComputeDRDYFunction(TS ts, PetscReal t, Vec U, Vec *DRDU)
1086d71ae5a4SJacob Faibussowitsch {
1087814e85d6SHong Zhang   PetscFunctionBegin;
1088814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1089c9ad9de0SHong Zhang   PetscValidHeaderSpecific(U, VEC_CLASSID, 3);
1090814e85d6SHong Zhang 
1091792fecdfSBarry Smith   PetscCallBack("TS callback DRDY for sensitivity analysis", (*ts->drdufunction)(ts, t, U, DRDU, ts->costintegrandctx));
10923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1093814e85d6SHong Zhang }
1094814e85d6SHong Zhang 
1095814e85d6SHong Zhang /*@
1096bcf0153eSBarry Smith   TSAdjointComputeDRDPFunction - Deprecated, use `TSGetQuadratureTS()` then `TSComputeRHSJacobianP()`
1097814e85d6SHong Zhang 
1098814e85d6SHong Zhang   Level: deprecated
1099814e85d6SHong Zhang 
1100814e85d6SHong Zhang @*/
1101d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointComputeDRDPFunction(TS ts, PetscReal t, Vec U, Vec *DRDP)
1102d71ae5a4SJacob Faibussowitsch {
1103814e85d6SHong Zhang   PetscFunctionBegin;
1104814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1105c9ad9de0SHong Zhang   PetscValidHeaderSpecific(U, VEC_CLASSID, 3);
1106814e85d6SHong Zhang 
1107792fecdfSBarry Smith   PetscCallBack("TS callback DRDP for sensitivity analysis", (*ts->drdpfunction)(ts, t, U, DRDP, ts->costintegrandctx));
11083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1109814e85d6SHong Zhang }
1110814e85d6SHong Zhang 
1111814e85d6SHong Zhang /*@C
1112814e85d6SHong Zhang    TSAdjointMonitorSensi - monitors the first lambda sensitivity
1113814e85d6SHong Zhang 
1114814e85d6SHong Zhang    Level: intermediate
1115814e85d6SHong Zhang 
1116bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSAdjointMonitorSet()`
1117814e85d6SHong Zhang @*/
1118d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointMonitorSensi(TS ts, PetscInt step, PetscReal ptime, Vec v, PetscInt numcost, Vec *lambda, Vec *mu, PetscViewerAndFormat *vf)
1119d71ae5a4SJacob Faibussowitsch {
1120814e85d6SHong Zhang   PetscViewer viewer = vf->viewer;
1121814e85d6SHong Zhang 
1122814e85d6SHong Zhang   PetscFunctionBegin;
1123064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 8);
11249566063dSJacob Faibussowitsch   PetscCall(PetscViewerPushFormat(viewer, vf->format));
11259566063dSJacob Faibussowitsch   PetscCall(VecView(lambda[0], viewer));
11269566063dSJacob Faibussowitsch   PetscCall(PetscViewerPopFormat(viewer));
11273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1128814e85d6SHong Zhang }
1129814e85d6SHong Zhang 
1130814e85d6SHong Zhang /*@C
1131814e85d6SHong Zhang    TSAdjointMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user
1132814e85d6SHong Zhang 
1133c3339decSBarry Smith    Collective
1134814e85d6SHong Zhang 
1135814e85d6SHong Zhang    Input Parameters:
1136bcf0153eSBarry Smith +  ts - `TS` object you wish to monitor
1137814e85d6SHong Zhang .  name - the monitor type one is seeking
1138814e85d6SHong Zhang .  help - message indicating what monitoring is done
1139814e85d6SHong Zhang .  manual - manual page for the monitor
1140814e85d6SHong Zhang .  monitor - the monitor function
1141bcf0153eSBarry Smith -  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
1142814e85d6SHong Zhang 
1143814e85d6SHong Zhang    Level: developer
1144814e85d6SHong Zhang 
1145bcf0153eSBarry Smith .seealso: [](chapter_ts), `PetscOptionsGetViewer()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`,
1146db781477SPatrick Sanan           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`
1147db781477SPatrick Sanan           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`,
1148db781477SPatrick Sanan           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
1149c2e3fba1SPatrick Sanan           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
1150db781477SPatrick Sanan           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
1151db781477SPatrick Sanan           `PetscOptionsFList()`, `PetscOptionsEList()`
1152814e85d6SHong Zhang @*/
1153d71ae5a4SJacob Faibussowitsch 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 *))
1154d71ae5a4SJacob Faibussowitsch {
1155814e85d6SHong Zhang   PetscViewer       viewer;
1156814e85d6SHong Zhang   PetscViewerFormat format;
1157814e85d6SHong Zhang   PetscBool         flg;
1158814e85d6SHong Zhang 
1159814e85d6SHong Zhang   PetscFunctionBegin;
11609566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)ts), ((PetscObject)ts)->options, ((PetscObject)ts)->prefix, name, &viewer, &format, &flg));
1161814e85d6SHong Zhang   if (flg) {
1162814e85d6SHong Zhang     PetscViewerAndFormat *vf;
11639566063dSJacob Faibussowitsch     PetscCall(PetscViewerAndFormatCreate(viewer, format, &vf));
11649566063dSJacob Faibussowitsch     PetscCall(PetscObjectDereference((PetscObject)viewer));
11651baa6e33SBarry Smith     if (monitorsetup) PetscCall((*monitorsetup)(ts, vf));
11669566063dSJacob Faibussowitsch     PetscCall(TSAdjointMonitorSet(ts, (PetscErrorCode(*)(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *, void *))monitor, vf, (PetscErrorCode(*)(void **))PetscViewerAndFormatDestroy));
1167814e85d6SHong Zhang   }
11683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1169814e85d6SHong Zhang }
1170814e85d6SHong Zhang 
1171814e85d6SHong Zhang /*@C
1172814e85d6SHong Zhang    TSAdjointMonitorSet - Sets an ADDITIONAL function that is to be used at every
1173814e85d6SHong Zhang    timestep to display the iteration's  progress.
1174814e85d6SHong Zhang 
1175c3339decSBarry Smith    Logically Collective
1176814e85d6SHong Zhang 
1177814e85d6SHong Zhang    Input Parameters:
1178bcf0153eSBarry Smith +  ts - the `TS` context obtained from `TSCreate()`
1179814e85d6SHong Zhang .  adjointmonitor - monitoring routine
1180814e85d6SHong Zhang .  adjointmctx - [optional] user-defined context for private data for the
1181814e85d6SHong Zhang              monitor routine (use NULL if no context is desired)
1182814e85d6SHong Zhang -  adjointmonitordestroy - [optional] routine that frees monitor context
1183814e85d6SHong Zhang           (may be NULL)
1184814e85d6SHong Zhang 
1185814e85d6SHong Zhang    Calling sequence of monitor:
1186814e85d6SHong Zhang $    int adjointmonitor(TS ts,PetscInt steps,PetscReal time,Vec u,PetscInt numcost,Vec *lambda, Vec *mu,void *adjointmctx)
1187814e85d6SHong Zhang 
1188bcf0153eSBarry Smith +    ts - the `TS` context
1189814e85d6SHong 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
1190814e85d6SHong Zhang                                been interpolated to)
1191814e85d6SHong Zhang .    time - current time
1192814e85d6SHong Zhang .    u - current iterate
1193814e85d6SHong Zhang .    numcost - number of cost functionos
1194814e85d6SHong Zhang .    lambda - sensitivities to initial conditions
1195814e85d6SHong Zhang .    mu - sensitivities to parameters
1196814e85d6SHong Zhang -    adjointmctx - [optional] adjoint monitoring context
1197814e85d6SHong Zhang 
1198bcf0153eSBarry Smith    Level: intermediate
1199bcf0153eSBarry Smith 
1200bcf0153eSBarry Smith    Note:
1201814e85d6SHong Zhang    This routine adds an additional monitor to the list of monitors that
1202814e85d6SHong Zhang    already has been loaded.
1203814e85d6SHong Zhang 
1204bcf0153eSBarry Smith    Fortran Note:
1205814e85d6SHong Zhang    Only a single monitor function can be set for each TS object
1206814e85d6SHong Zhang 
1207bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSAdjointSolve()`, `TSAdjointMonitorCancel()`
1208814e85d6SHong Zhang @*/
1209d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointMonitorSet(TS ts, PetscErrorCode (*adjointmonitor)(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *, void *), void *adjointmctx, PetscErrorCode (*adjointmdestroy)(void **))
1210d71ae5a4SJacob Faibussowitsch {
1211814e85d6SHong Zhang   PetscInt  i;
1212814e85d6SHong Zhang   PetscBool identical;
1213814e85d6SHong Zhang 
1214814e85d6SHong Zhang   PetscFunctionBegin;
1215814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1216814e85d6SHong Zhang   for (i = 0; i < ts->numbermonitors; i++) {
12179566063dSJacob Faibussowitsch     PetscCall(PetscMonitorCompare((PetscErrorCode(*)(void))adjointmonitor, adjointmctx, adjointmdestroy, (PetscErrorCode(*)(void))ts->adjointmonitor[i], ts->adjointmonitorcontext[i], ts->adjointmonitordestroy[i], &identical));
12183ba16761SJacob Faibussowitsch     if (identical) PetscFunctionReturn(PETSC_SUCCESS);
1219814e85d6SHong Zhang   }
12203c633725SBarry Smith   PetscCheck(ts->numberadjointmonitors < MAXTSMONITORS, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many adjoint monitors set");
1221814e85d6SHong Zhang   ts->adjointmonitor[ts->numberadjointmonitors]          = adjointmonitor;
1222814e85d6SHong Zhang   ts->adjointmonitordestroy[ts->numberadjointmonitors]   = adjointmdestroy;
1223814e85d6SHong Zhang   ts->adjointmonitorcontext[ts->numberadjointmonitors++] = (void *)adjointmctx;
12243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1225814e85d6SHong Zhang }
1226814e85d6SHong Zhang 
1227814e85d6SHong Zhang /*@C
1228814e85d6SHong Zhang    TSAdjointMonitorCancel - Clears all the adjoint monitors that have been set on a time-step object.
1229814e85d6SHong Zhang 
1230c3339decSBarry Smith    Logically Collective
1231814e85d6SHong Zhang 
1232814e85d6SHong Zhang    Input Parameters:
1233bcf0153eSBarry Smith .  ts - the `TS` context obtained from `TSCreate()`
1234814e85d6SHong Zhang 
1235814e85d6SHong Zhang    Notes:
1236814e85d6SHong Zhang    There is no way to remove a single, specific monitor.
1237814e85d6SHong Zhang 
1238814e85d6SHong Zhang    Level: intermediate
1239814e85d6SHong Zhang 
1240bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSAdjointSolve()`, `TSAdjointMonitorSet()`
1241814e85d6SHong Zhang @*/
1242d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointMonitorCancel(TS ts)
1243d71ae5a4SJacob Faibussowitsch {
1244814e85d6SHong Zhang   PetscInt i;
1245814e85d6SHong Zhang 
1246814e85d6SHong Zhang   PetscFunctionBegin;
1247814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1248814e85d6SHong Zhang   for (i = 0; i < ts->numberadjointmonitors; i++) {
124948a46eb9SPierre Jolivet     if (ts->adjointmonitordestroy[i]) PetscCall((*ts->adjointmonitordestroy[i])(&ts->adjointmonitorcontext[i]));
1250814e85d6SHong Zhang   }
1251814e85d6SHong Zhang   ts->numberadjointmonitors = 0;
12523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1253814e85d6SHong Zhang }
1254814e85d6SHong Zhang 
1255814e85d6SHong Zhang /*@C
1256814e85d6SHong Zhang    TSAdjointMonitorDefault - the default monitor of adjoint computations
1257814e85d6SHong Zhang 
1258814e85d6SHong Zhang    Level: intermediate
1259814e85d6SHong Zhang 
1260bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSAdjointSolve()`, `TSAdjointMonitorSet()`
1261814e85d6SHong Zhang @*/
1262d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointMonitorDefault(TS ts, PetscInt step, PetscReal ptime, Vec v, PetscInt numcost, Vec *lambda, Vec *mu, PetscViewerAndFormat *vf)
1263d71ae5a4SJacob Faibussowitsch {
1264814e85d6SHong Zhang   PetscViewer viewer = vf->viewer;
1265814e85d6SHong Zhang 
1266814e85d6SHong Zhang   PetscFunctionBegin;
1267064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 8);
12689566063dSJacob Faibussowitsch   PetscCall(PetscViewerPushFormat(viewer, vf->format));
12699566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)ts)->tablevel));
127063a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " TS dt %g time %g%s", step, (double)ts->time_step, (double)ptime, ts->steprollback ? " (r)\n" : "\n"));
12719566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)ts)->tablevel));
12729566063dSJacob Faibussowitsch   PetscCall(PetscViewerPopFormat(viewer));
12733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1274814e85d6SHong Zhang }
1275814e85d6SHong Zhang 
1276814e85d6SHong Zhang /*@C
1277bcf0153eSBarry Smith    TSAdjointMonitorDrawSensi - Monitors progress of the adjoint `TS` solvers by calling
1278bcf0153eSBarry Smith    `VecView()` for the sensitivities to initial states at each timestep
1279814e85d6SHong Zhang 
1280c3339decSBarry Smith    Collective
1281814e85d6SHong Zhang 
1282814e85d6SHong Zhang    Input Parameters:
1283bcf0153eSBarry Smith +  ts - the `TS` context
1284814e85d6SHong Zhang .  step - current time-step
1285814e85d6SHong Zhang .  ptime - current time
1286814e85d6SHong Zhang .  u - current state
1287814e85d6SHong Zhang .  numcost - number of cost functions
1288814e85d6SHong Zhang .  lambda - sensitivities to initial conditions
1289814e85d6SHong Zhang .  mu - sensitivities to parameters
1290814e85d6SHong Zhang -  dummy - either a viewer or NULL
1291814e85d6SHong Zhang 
1292814e85d6SHong Zhang    Level: intermediate
1293814e85d6SHong Zhang 
1294bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSAdjointSolve()`, `TSAdjointMonitorSet()`, `TSAdjointMonitorDefault()`, `VecView()`
1295814e85d6SHong Zhang @*/
1296d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointMonitorDrawSensi(TS ts, PetscInt step, PetscReal ptime, Vec u, PetscInt numcost, Vec *lambda, Vec *mu, void *dummy)
1297d71ae5a4SJacob Faibussowitsch {
1298814e85d6SHong Zhang   TSMonitorDrawCtx ictx = (TSMonitorDrawCtx)dummy;
1299814e85d6SHong Zhang   PetscDraw        draw;
1300814e85d6SHong Zhang   PetscReal        xl, yl, xr, yr, h;
1301814e85d6SHong Zhang   char             time[32];
1302814e85d6SHong Zhang 
1303814e85d6SHong Zhang   PetscFunctionBegin;
13043ba16761SJacob Faibussowitsch   if (!(((ictx->howoften > 0) && (!(step % ictx->howoften))) || ((ictx->howoften == -1) && ts->reason))) PetscFunctionReturn(PETSC_SUCCESS);
1305814e85d6SHong Zhang 
13069566063dSJacob Faibussowitsch   PetscCall(VecView(lambda[0], ictx->viewer));
13079566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawGetDraw(ictx->viewer, 0, &draw));
13089566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(time, 32, "Timestep %d Time %g", (int)step, (double)ptime));
13099566063dSJacob Faibussowitsch   PetscCall(PetscDrawGetCoordinates(draw, &xl, &yl, &xr, &yr));
1310814e85d6SHong Zhang   h = yl + .95 * (yr - yl);
13119566063dSJacob Faibussowitsch   PetscCall(PetscDrawStringCentered(draw, .5 * (xl + xr), h, PETSC_DRAW_BLACK, time));
13129566063dSJacob Faibussowitsch   PetscCall(PetscDrawFlush(draw));
13133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1314814e85d6SHong Zhang }
1315814e85d6SHong Zhang 
1316814e85d6SHong Zhang /*
1317bcf0153eSBarry Smith    TSAdjointSetFromOptions - Sets various `TS` adjoint parameters from user options.
1318814e85d6SHong Zhang 
1319c3339decSBarry Smith    Collective
1320814e85d6SHong Zhang 
1321814e85d6SHong Zhang    Input Parameter:
1322bcf0153eSBarry Smith .  ts - the `TS` context
1323814e85d6SHong Zhang 
1324814e85d6SHong Zhang    Options Database Keys:
1325814e85d6SHong Zhang +  -ts_adjoint_solve <yes,no> After solving the ODE/DAE solve the adjoint problem (requires -ts_save_trajectory)
1326814e85d6SHong Zhang .  -ts_adjoint_monitor - print information at each adjoint time step
1327814e85d6SHong Zhang -  -ts_adjoint_monitor_draw_sensi - monitor the sensitivity of the first cost function wrt initial conditions (lambda[0]) graphically
1328814e85d6SHong Zhang 
1329814e85d6SHong Zhang    Level: developer
1330814e85d6SHong Zhang 
1331bcf0153eSBarry Smith    Note:
1332814e85d6SHong Zhang     This is not normally called directly by users
1333814e85d6SHong Zhang 
1334bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSSetSaveTrajectory()`, `TSTrajectorySetUp()`
1335814e85d6SHong Zhang */
1336d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointSetFromOptions(TS ts, PetscOptionItems *PetscOptionsObject)
1337d71ae5a4SJacob Faibussowitsch {
1338814e85d6SHong Zhang   PetscBool tflg, opt;
1339814e85d6SHong Zhang 
1340814e85d6SHong Zhang   PetscFunctionBegin;
1341dbbe0bcdSBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1342d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "TS Adjoint options");
1343814e85d6SHong Zhang   tflg = ts->adjoint_solve ? PETSC_TRUE : PETSC_FALSE;
13449566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-ts_adjoint_solve", "Solve the adjoint problem immediately after solving the forward problem", "", tflg, &tflg, &opt));
1345814e85d6SHong Zhang   if (opt) {
13469566063dSJacob Faibussowitsch     PetscCall(TSSetSaveTrajectory(ts));
1347814e85d6SHong Zhang     ts->adjoint_solve = tflg;
1348814e85d6SHong Zhang   }
13499566063dSJacob Faibussowitsch   PetscCall(TSAdjointMonitorSetFromOptions(ts, "-ts_adjoint_monitor", "Monitor adjoint timestep size", "TSAdjointMonitorDefault", TSAdjointMonitorDefault, NULL));
13509566063dSJacob Faibussowitsch   PetscCall(TSAdjointMonitorSetFromOptions(ts, "-ts_adjoint_monitor_sensi", "Monitor sensitivity in the adjoint computation", "TSAdjointMonitorSensi", TSAdjointMonitorSensi, NULL));
1351814e85d6SHong Zhang   opt = PETSC_FALSE;
13529566063dSJacob Faibussowitsch   PetscCall(PetscOptionsName("-ts_adjoint_monitor_draw_sensi", "Monitor adjoint sensitivities (lambda only) graphically", "TSAdjointMonitorDrawSensi", &opt));
1353814e85d6SHong Zhang   if (opt) {
1354814e85d6SHong Zhang     TSMonitorDrawCtx ctx;
1355814e85d6SHong Zhang     PetscInt         howoften = 1;
1356814e85d6SHong Zhang 
13579566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-ts_adjoint_monitor_draw_sensi", "Monitor adjoint sensitivities (lambda only) graphically", "TSAdjointMonitorDrawSensi", howoften, &howoften, NULL));
13589566063dSJacob Faibussowitsch     PetscCall(TSMonitorDrawCtxCreate(PetscObjectComm((PetscObject)ts), NULL, NULL, PETSC_DECIDE, PETSC_DECIDE, 300, 300, howoften, &ctx));
13599566063dSJacob Faibussowitsch     PetscCall(TSAdjointMonitorSet(ts, TSAdjointMonitorDrawSensi, ctx, (PetscErrorCode(*)(void **))TSMonitorDrawCtxDestroy));
1360814e85d6SHong Zhang   }
13613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1362814e85d6SHong Zhang }
1363814e85d6SHong Zhang 
1364814e85d6SHong Zhang /*@
1365814e85d6SHong Zhang    TSAdjointStep - Steps one time step backward in the adjoint run
1366814e85d6SHong Zhang 
1367c3339decSBarry Smith    Collective
1368814e85d6SHong Zhang 
1369814e85d6SHong Zhang    Input Parameter:
1370bcf0153eSBarry Smith .  ts - the `TS` context obtained from `TSCreate()`
1371814e85d6SHong Zhang 
1372814e85d6SHong Zhang    Level: intermediate
1373814e85d6SHong Zhang 
1374bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSAdjointSetUp()`, `TSAdjointSolve()`
1375814e85d6SHong Zhang @*/
1376d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointStep(TS ts)
1377d71ae5a4SJacob Faibussowitsch {
1378814e85d6SHong Zhang   DM dm;
1379814e85d6SHong Zhang 
1380814e85d6SHong Zhang   PetscFunctionBegin;
1381814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
13829566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
13839566063dSJacob Faibussowitsch   PetscCall(TSAdjointSetUp(ts));
13847b0e2f17SHong Zhang   ts->steps--; /* must decrease the step index before the adjoint step is taken. */
1385814e85d6SHong Zhang 
1386814e85d6SHong Zhang   ts->reason     = TS_CONVERGED_ITERATING;
1387814e85d6SHong Zhang   ts->ptime_prev = ts->ptime;
13889566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(TS_AdjointStep, ts, 0, 0, 0));
1389dbbe0bcdSBarry Smith   PetscUseTypeMethod(ts, adjointstep);
13909566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(TS_AdjointStep, ts, 0, 0, 0));
13917b0e2f17SHong Zhang   ts->adjoint_steps++;
1392814e85d6SHong Zhang 
1393814e85d6SHong Zhang   if (ts->reason < 0) {
13943c633725SBarry Smith     PetscCheck(!ts->errorifstepfailed, PetscObjectComm((PetscObject)ts), PETSC_ERR_NOT_CONVERGED, "TSAdjointStep has failed due to %s", TSConvergedReasons[ts->reason]);
1395814e85d6SHong Zhang   } else if (!ts->reason) {
1396814e85d6SHong Zhang     if (ts->adjoint_steps >= ts->adjoint_max_steps) ts->reason = TS_CONVERGED_ITS;
1397814e85d6SHong Zhang   }
13983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1399814e85d6SHong Zhang }
1400814e85d6SHong Zhang 
1401814e85d6SHong Zhang /*@
1402814e85d6SHong Zhang    TSAdjointSolve - Solves the discrete ajoint problem for an ODE/DAE
1403814e85d6SHong Zhang 
1404c3339decSBarry Smith    Collective
1405bcf0153eSBarry Smith `
1406814e85d6SHong Zhang    Input Parameter:
1407bcf0153eSBarry Smith .  ts - the `TS` context obtained from `TSCreate()`
1408814e85d6SHong Zhang 
1409bcf0153eSBarry Smith    Options Database Key:
1410814e85d6SHong Zhang . -ts_adjoint_view_solution <viewerinfo> - views the first gradient with respect to the initial values
1411814e85d6SHong Zhang 
1412814e85d6SHong Zhang    Level: intermediate
1413814e85d6SHong Zhang 
1414814e85d6SHong Zhang    Notes:
1415bcf0153eSBarry Smith    This must be called after a call to `TSSolve()` that solves the forward problem
1416814e85d6SHong Zhang 
1417bcf0153eSBarry Smith    By default this will integrate back to the initial time, one can use `TSAdjointSetSteps()` to step back to a later time
1418814e85d6SHong Zhang 
1419bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSAdjointSolve()`, `TSCreate()`, `TSSetCostGradients()`, `TSSetSolution()`, `TSAdjointStep()`
1420814e85d6SHong Zhang @*/
1421d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointSolve(TS ts)
1422d71ae5a4SJacob Faibussowitsch {
1423f1d62c27SHong Zhang   static PetscBool cite = PETSC_FALSE;
14247f59fb53SHong Zhang #if defined(TSADJOINT_STAGE)
14257f59fb53SHong Zhang   PetscLogStage adjoint_stage;
14267f59fb53SHong Zhang #endif
1427814e85d6SHong Zhang 
1428814e85d6SHong Zhang   PetscFunctionBegin;
1429814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1430421b783aSHong Zhang   PetscCall(PetscCitationsRegister("@article{Zhang2022tsadjoint,\n"
1431421b783aSHong Zhang                                    "  title         = {{PETSc TSAdjoint: A Discrete Adjoint ODE Solver for First-Order and Second-Order Sensitivity Analysis}},\n"
1432f1d62c27SHong Zhang                                    "  author        = {Zhang, Hong and Constantinescu, Emil M.  and Smith, Barry F.},\n"
1433421b783aSHong Zhang                                    "  journal       = {SIAM Journal on Scientific Computing},\n"
1434421b783aSHong Zhang                                    "  volume        = {44},\n"
1435421b783aSHong Zhang                                    "  number        = {1},\n"
1436421b783aSHong Zhang                                    "  pages         = {C1-C24},\n"
1437421b783aSHong Zhang                                    "  doi           = {10.1137/21M140078X},\n"
14389371c9d4SSatish Balay                                    "  year          = {2022}\n}\n",
14399371c9d4SSatish Balay                                    &cite));
14407f59fb53SHong Zhang #if defined(TSADJOINT_STAGE)
14419566063dSJacob Faibussowitsch   PetscCall(PetscLogStageRegister("TSAdjoint", &adjoint_stage));
14429566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePush(adjoint_stage));
14437f59fb53SHong Zhang #endif
14449566063dSJacob Faibussowitsch   PetscCall(TSAdjointSetUp(ts));
1445814e85d6SHong Zhang 
1446814e85d6SHong Zhang   /* reset time step and iteration counters */
1447814e85d6SHong Zhang   ts->adjoint_steps     = 0;
1448814e85d6SHong Zhang   ts->ksp_its           = 0;
1449814e85d6SHong Zhang   ts->snes_its          = 0;
1450814e85d6SHong Zhang   ts->num_snes_failures = 0;
1451814e85d6SHong Zhang   ts->reject            = 0;
1452814e85d6SHong Zhang   ts->reason            = TS_CONVERGED_ITERATING;
1453814e85d6SHong Zhang 
1454814e85d6SHong Zhang   if (!ts->adjoint_max_steps) ts->adjoint_max_steps = ts->steps;
1455814e85d6SHong Zhang   if (ts->adjoint_steps >= ts->adjoint_max_steps) ts->reason = TS_CONVERGED_ITS;
1456814e85d6SHong Zhang 
1457814e85d6SHong Zhang   while (!ts->reason) {
14589566063dSJacob Faibussowitsch     PetscCall(TSTrajectoryGet(ts->trajectory, ts, ts->steps, &ts->ptime));
14599566063dSJacob Faibussowitsch     PetscCall(TSAdjointMonitor(ts, ts->steps, ts->ptime, ts->vec_sol, ts->numcost, ts->vecs_sensi, ts->vecs_sensip));
14609566063dSJacob Faibussowitsch     PetscCall(TSAdjointEventHandler(ts));
14619566063dSJacob Faibussowitsch     PetscCall(TSAdjointStep(ts));
146248a46eb9SPierre Jolivet     if ((ts->vec_costintegral || ts->quadraturets) && !ts->costintegralfwd) PetscCall(TSAdjointCostIntegral(ts));
1463814e85d6SHong Zhang   }
1464bdbff044SHong Zhang   if (!ts->steps) {
14659566063dSJacob Faibussowitsch     PetscCall(TSTrajectoryGet(ts->trajectory, ts, ts->steps, &ts->ptime));
14669566063dSJacob Faibussowitsch     PetscCall(TSAdjointMonitor(ts, ts->steps, ts->ptime, ts->vec_sol, ts->numcost, ts->vecs_sensi, ts->vecs_sensip));
1467bdbff044SHong Zhang   }
1468814e85d6SHong Zhang   ts->solvetime = ts->ptime;
14699566063dSJacob Faibussowitsch   PetscCall(TSTrajectoryViewFromOptions(ts->trajectory, NULL, "-ts_trajectory_view"));
14709566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(ts->vecs_sensi[0], (PetscObject)ts, "-ts_adjoint_view_solution"));
1471814e85d6SHong Zhang   ts->adjoint_max_steps = 0;
14727f59fb53SHong Zhang #if defined(TSADJOINT_STAGE)
14739566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePop());
14747f59fb53SHong Zhang #endif
14753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1476814e85d6SHong Zhang }
1477814e85d6SHong Zhang 
1478814e85d6SHong Zhang /*@C
1479bcf0153eSBarry Smith    TSAdjointMonitor - Runs all user-provided adjoint monitor routines set using `TSAdjointMonitorSet()`
1480814e85d6SHong Zhang 
1481c3339decSBarry Smith    Collective
1482814e85d6SHong Zhang 
1483814e85d6SHong Zhang    Input Parameters:
1484bcf0153eSBarry Smith +  ts - time stepping context obtained from `TSCreate()`
1485814e85d6SHong Zhang .  step - step number that has just completed
1486814e85d6SHong Zhang .  ptime - model time of the state
1487814e85d6SHong Zhang .  u - state at the current model time
1488814e85d6SHong Zhang .  numcost - number of cost functions (dimension of lambda  or mu)
1489814e85d6SHong Zhang .  lambda - vectors containing the gradients of the cost functions with respect to the ODE/DAE solution variables
1490814e85d6SHong Zhang -  mu - vectors containing the gradients of the cost functions with respect to the problem parameters
1491814e85d6SHong Zhang 
1492814e85d6SHong Zhang    Level: developer
1493814e85d6SHong Zhang 
1494bcf0153eSBarry Smith    Note:
1495bcf0153eSBarry Smith    `TSAdjointMonitor()` is typically used automatically within the time stepping implementations.
1496bcf0153eSBarry Smith    Users would almost never call this routine directly.
1497bcf0153eSBarry Smith 
1498bcf0153eSBarry Smith .seealso: `TSAdjointMonitorSet()`, `TSAdjointSolve()`
1499814e85d6SHong Zhang @*/
1500d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointMonitor(TS ts, PetscInt step, PetscReal ptime, Vec u, PetscInt numcost, Vec *lambda, Vec *mu)
1501d71ae5a4SJacob Faibussowitsch {
1502814e85d6SHong Zhang   PetscInt i, n = ts->numberadjointmonitors;
1503814e85d6SHong Zhang 
1504814e85d6SHong Zhang   PetscFunctionBegin;
1505814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1506814e85d6SHong Zhang   PetscValidHeaderSpecific(u, VEC_CLASSID, 4);
15079566063dSJacob Faibussowitsch   PetscCall(VecLockReadPush(u));
150848a46eb9SPierre Jolivet   for (i = 0; i < n; i++) PetscCall((*ts->adjointmonitor[i])(ts, step, ptime, u, numcost, lambda, mu, ts->adjointmonitorcontext[i]));
15099566063dSJacob Faibussowitsch   PetscCall(VecLockReadPop(u));
15103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1511814e85d6SHong Zhang }
1512814e85d6SHong Zhang 
1513814e85d6SHong Zhang /*@
1514814e85d6SHong Zhang  TSAdjointCostIntegral - Evaluate the cost integral in the adjoint run.
1515814e85d6SHong Zhang 
1516c3339decSBarry Smith  Collective
1517814e85d6SHong Zhang 
15184165533cSJose E. Roman  Input Parameter:
1519814e85d6SHong Zhang  .  ts - time stepping context
1520814e85d6SHong Zhang 
1521814e85d6SHong Zhang  Level: advanced
1522814e85d6SHong Zhang 
1523814e85d6SHong Zhang  Notes:
1524bcf0153eSBarry Smith  This function cannot be called until `TSAdjointStep()` has been completed.
1525814e85d6SHong Zhang 
1526bcf0153eSBarry Smith  .seealso: [](chapter_ts), `TSAdjointSolve()`, `TSAdjointStep()`
1527814e85d6SHong Zhang  @*/
1528d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointCostIntegral(TS ts)
1529d71ae5a4SJacob Faibussowitsch {
1530362febeeSStefano Zampini   PetscFunctionBegin;
1531814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1532dbbe0bcdSBarry Smith   PetscUseTypeMethod(ts, adjointintegral);
15333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1534814e85d6SHong Zhang }
1535814e85d6SHong Zhang 
1536814e85d6SHong Zhang /* ------------------ Forward (tangent linear) sensitivity  ------------------*/
1537814e85d6SHong Zhang 
1538814e85d6SHong Zhang /*@
1539814e85d6SHong Zhang   TSForwardSetUp - Sets up the internal data structures for the later use
1540814e85d6SHong Zhang   of forward sensitivity analysis
1541814e85d6SHong Zhang 
1542c3339decSBarry Smith   Collective
1543814e85d6SHong Zhang 
1544814e85d6SHong Zhang   Input Parameter:
1545bcf0153eSBarry Smith . ts - the `TS` context obtained from `TSCreate()`
1546814e85d6SHong Zhang 
1547814e85d6SHong Zhang   Level: advanced
1548814e85d6SHong Zhang 
1549bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSCreate()`, `TSDestroy()`, `TSSetUp()`
1550814e85d6SHong Zhang @*/
1551d71ae5a4SJacob Faibussowitsch PetscErrorCode TSForwardSetUp(TS ts)
1552d71ae5a4SJacob Faibussowitsch {
1553814e85d6SHong Zhang   PetscFunctionBegin;
1554814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
15553ba16761SJacob Faibussowitsch   if (ts->forwardsetupcalled) PetscFunctionReturn(PETSC_SUCCESS);
1556dbbe0bcdSBarry Smith   PetscTryTypeMethod(ts, forwardsetup);
15579566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &ts->vec_sensip_col));
1558814e85d6SHong Zhang   ts->forwardsetupcalled = PETSC_TRUE;
15593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1560814e85d6SHong Zhang }
1561814e85d6SHong Zhang 
1562814e85d6SHong Zhang /*@
15637adebddeSHong Zhang   TSForwardReset - Reset the internal data structures used by forward sensitivity analysis
15647adebddeSHong Zhang 
1565c3339decSBarry Smith   Collective
15667adebddeSHong Zhang 
15677adebddeSHong Zhang   Input Parameter:
1568bcf0153eSBarry Smith . ts - the `TS` context obtained from `TSCreate()`
15697adebddeSHong Zhang 
15707adebddeSHong Zhang   Level: advanced
15717adebddeSHong Zhang 
1572bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSCreate()`, `TSDestroy()`, `TSForwardSetUp()`
15737adebddeSHong Zhang @*/
1574d71ae5a4SJacob Faibussowitsch PetscErrorCode TSForwardReset(TS ts)
1575d71ae5a4SJacob Faibussowitsch {
15767207e0fdSHong Zhang   TS quadts = ts->quadraturets;
15777adebddeSHong Zhang 
15787adebddeSHong Zhang   PetscFunctionBegin;
15797adebddeSHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1580dbbe0bcdSBarry Smith   PetscTryTypeMethod(ts, forwardreset);
15819566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&ts->mat_sensip));
158248a46eb9SPierre Jolivet   if (quadts) PetscCall(MatDestroy(&quadts->mat_sensip));
15839566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ts->vec_sensip_col));
15847adebddeSHong Zhang   ts->forward_solve      = PETSC_FALSE;
15857adebddeSHong Zhang   ts->forwardsetupcalled = PETSC_FALSE;
15863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15877adebddeSHong Zhang }
15887adebddeSHong Zhang 
15897adebddeSHong Zhang /*@
1590814e85d6SHong Zhang   TSForwardSetIntegralGradients - Set the vectors holding forward sensitivities of the integral term.
1591814e85d6SHong Zhang 
1592d8d19677SJose E. Roman   Input Parameters:
1593bcf0153eSBarry Smith + ts - the `TS` context obtained from `TSCreate()`
1594814e85d6SHong Zhang . numfwdint - number of integrals
159567b8a455SSatish Balay - vp - the vectors containing the gradients for each integral w.r.t. parameters
1596814e85d6SHong Zhang 
15977207e0fdSHong Zhang   Level: deprecated
1598814e85d6SHong Zhang 
1599bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSForwardGetSensitivities()`, `TSForwardSetIntegralGradients()`, `TSForwardGetIntegralGradients()`, `TSForwardStep()`
1600814e85d6SHong Zhang @*/
1601d71ae5a4SJacob Faibussowitsch PetscErrorCode TSForwardSetIntegralGradients(TS ts, PetscInt numfwdint, Vec *vp)
1602d71ae5a4SJacob Faibussowitsch {
1603814e85d6SHong Zhang   PetscFunctionBegin;
1604814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
16053c633725SBarry Smith   PetscCheck(!ts->numcost || ts->numcost == numfwdint, PetscObjectComm((PetscObject)ts), PETSC_ERR_USER, "The number of cost functions (2nd parameter of TSSetCostIntegrand()) is inconsistent with the one set by TSSetCostIntegrand()");
1606814e85d6SHong Zhang   if (!ts->numcost) ts->numcost = numfwdint;
1607814e85d6SHong Zhang 
1608814e85d6SHong Zhang   ts->vecs_integral_sensip = vp;
16093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1610814e85d6SHong Zhang }
1611814e85d6SHong Zhang 
1612814e85d6SHong Zhang /*@
1613814e85d6SHong Zhang   TSForwardGetIntegralGradients - Returns the forward sensitivities ofthe integral term.
1614814e85d6SHong Zhang 
1615814e85d6SHong Zhang   Input Parameter:
1616bcf0153eSBarry Smith . ts - the `TS` context obtained from `TSCreate()`
1617814e85d6SHong Zhang 
1618814e85d6SHong Zhang   Output Parameter:
161967b8a455SSatish Balay . vp - the vectors containing the gradients for each integral w.r.t. parameters
1620814e85d6SHong Zhang 
16217207e0fdSHong Zhang   Level: deprecated
1622814e85d6SHong Zhang 
1623bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSForwardSetSensitivities()`, `TSForwardSetIntegralGradients()`, `TSForwardGetIntegralGradients()`, `TSForwardStep()`
1624814e85d6SHong Zhang @*/
1625d71ae5a4SJacob Faibussowitsch PetscErrorCode TSForwardGetIntegralGradients(TS ts, PetscInt *numfwdint, Vec **vp)
1626d71ae5a4SJacob Faibussowitsch {
1627814e85d6SHong Zhang   PetscFunctionBegin;
1628814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1629814e85d6SHong Zhang   PetscValidPointer(vp, 3);
1630814e85d6SHong Zhang   if (numfwdint) *numfwdint = ts->numcost;
1631814e85d6SHong Zhang   if (vp) *vp = ts->vecs_integral_sensip;
16323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1633814e85d6SHong Zhang }
1634814e85d6SHong Zhang 
1635814e85d6SHong Zhang /*@
1636814e85d6SHong Zhang   TSForwardStep - Compute the forward sensitivity for one time step.
1637814e85d6SHong Zhang 
1638c3339decSBarry Smith   Collective
1639814e85d6SHong Zhang 
16404165533cSJose E. Roman   Input Parameter:
1641814e85d6SHong Zhang . ts - time stepping context
1642814e85d6SHong Zhang 
1643814e85d6SHong Zhang   Level: advanced
1644814e85d6SHong Zhang 
1645814e85d6SHong Zhang   Notes:
1646bcf0153eSBarry Smith   This function cannot be called until `TSStep()` has been completed.
1647814e85d6SHong Zhang 
1648bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSForwardSetSensitivities()`, `TSForwardGetSensitivities()`, `TSForwardSetIntegralGradients()`, `TSForwardGetIntegralGradients()`, `TSForwardSetUp()`
1649814e85d6SHong Zhang @*/
1650d71ae5a4SJacob Faibussowitsch PetscErrorCode TSForwardStep(TS ts)
1651d71ae5a4SJacob Faibussowitsch {
1652362febeeSStefano Zampini   PetscFunctionBegin;
1653814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
16549566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(TS_ForwardStep, ts, 0, 0, 0));
1655dbbe0bcdSBarry Smith   PetscUseTypeMethod(ts, forwardstep);
16569566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(TS_ForwardStep, ts, 0, 0, 0));
16573c633725SBarry Smith   PetscCheck(ts->reason >= 0 || !ts->errorifstepfailed, PetscObjectComm((PetscObject)ts), PETSC_ERR_NOT_CONVERGED, "TSFowardStep has failed due to %s", TSConvergedReasons[ts->reason]);
16583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1659814e85d6SHong Zhang }
1660814e85d6SHong Zhang 
1661814e85d6SHong Zhang /*@
1662814e85d6SHong Zhang   TSForwardSetSensitivities - Sets the initial value of the trajectory sensitivities of solution  w.r.t. the problem parameters and initial values.
1663814e85d6SHong Zhang 
1664c3339decSBarry Smith   Logically Collective
1665814e85d6SHong Zhang 
1666814e85d6SHong Zhang   Input Parameters:
1667bcf0153eSBarry Smith + ts - the `TS` context obtained from `TSCreate()`
1668814e85d6SHong Zhang . nump - number of parameters
1669814e85d6SHong Zhang - Smat - sensitivities with respect to the parameters, the number of entries in these vectors is the same as the number of parameters
1670814e85d6SHong Zhang 
1671814e85d6SHong Zhang   Level: beginner
1672814e85d6SHong Zhang 
1673814e85d6SHong Zhang   Notes:
1674814e85d6SHong Zhang   Forward sensitivity is also called 'trajectory sensitivity' in some fields such as power systems.
1675bcf0153eSBarry Smith   This function turns on a flag to trigger `TSSolve()` to compute forward sensitivities automatically.
1676bcf0153eSBarry Smith   You must call this function before `TSSolve()`.
1677814e85d6SHong Zhang   The entries in the sensitivity matrix must be correctly initialized with the values S = dy/dp|startingtime.
1678814e85d6SHong Zhang 
1679bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSForwardGetSensitivities()`, `TSForwardSetIntegralGradients()`, `TSForwardGetIntegralGradients()`, `TSForwardStep()`
1680814e85d6SHong Zhang @*/
1681d71ae5a4SJacob Faibussowitsch PetscErrorCode TSForwardSetSensitivities(TS ts, PetscInt nump, Mat Smat)
1682d71ae5a4SJacob Faibussowitsch {
1683814e85d6SHong Zhang   PetscFunctionBegin;
1684814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1685814e85d6SHong Zhang   PetscValidHeaderSpecific(Smat, MAT_CLASSID, 3);
1686814e85d6SHong Zhang   ts->forward_solve = PETSC_TRUE;
1687b5b4017aSHong Zhang   if (nump == PETSC_DEFAULT) {
16889566063dSJacob Faibussowitsch     PetscCall(MatGetSize(Smat, NULL, &ts->num_parameters));
1689b5b4017aSHong Zhang   } else ts->num_parameters = nump;
16909566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)Smat));
16919566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&ts->mat_sensip));
1692814e85d6SHong Zhang   ts->mat_sensip = Smat;
16933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1694814e85d6SHong Zhang }
1695814e85d6SHong Zhang 
1696814e85d6SHong Zhang /*@
1697814e85d6SHong Zhang   TSForwardGetSensitivities - Returns the trajectory sensitivities
1698814e85d6SHong Zhang 
1699bcf0153eSBarry Smith   Not Collective, but Smat returned is parallel if ts is parallel
1700814e85d6SHong Zhang 
1701d8d19677SJose E. Roman   Output Parameters:
1702bcf0153eSBarry Smith + ts - the `TS` context obtained from `TSCreate()`
1703814e85d6SHong Zhang . nump - number of parameters
1704814e85d6SHong Zhang - Smat - sensitivities with respect to the parameters, the number of entries in these vectors is the same as the number of parameters
1705814e85d6SHong Zhang 
1706814e85d6SHong Zhang   Level: intermediate
1707814e85d6SHong Zhang 
1708bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSForwardSetSensitivities()`, `TSForwardSetIntegralGradients()`, `TSForwardGetIntegralGradients()`, `TSForwardStep()`
1709814e85d6SHong Zhang @*/
1710d71ae5a4SJacob Faibussowitsch PetscErrorCode TSForwardGetSensitivities(TS ts, PetscInt *nump, Mat *Smat)
1711d71ae5a4SJacob Faibussowitsch {
1712814e85d6SHong Zhang   PetscFunctionBegin;
1713814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1714814e85d6SHong Zhang   if (nump) *nump = ts->num_parameters;
1715814e85d6SHong Zhang   if (Smat) *Smat = ts->mat_sensip;
17163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1717814e85d6SHong Zhang }
1718814e85d6SHong Zhang 
1719814e85d6SHong Zhang /*@
1720814e85d6SHong Zhang    TSForwardCostIntegral - Evaluate the cost integral in the forward run.
1721814e85d6SHong Zhang 
1722c3339decSBarry Smith    Collective
1723814e85d6SHong Zhang 
17244165533cSJose E. Roman    Input Parameter:
1725814e85d6SHong Zhang .  ts - time stepping context
1726814e85d6SHong Zhang 
1727814e85d6SHong Zhang    Level: advanced
1728814e85d6SHong Zhang 
1729bcf0153eSBarry Smith    Note:
1730bcf0153eSBarry Smith    This function cannot be called until `TSStep()` has been completed.
1731814e85d6SHong Zhang 
1732bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSSolve()`, `TSAdjointCostIntegral()`
1733814e85d6SHong Zhang @*/
1734d71ae5a4SJacob Faibussowitsch PetscErrorCode TSForwardCostIntegral(TS ts)
1735d71ae5a4SJacob Faibussowitsch {
1736362febeeSStefano Zampini   PetscFunctionBegin;
1737814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1738dbbe0bcdSBarry Smith   PetscUseTypeMethod(ts, forwardintegral);
17393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1740814e85d6SHong Zhang }
1741b5b4017aSHong Zhang 
1742b5b4017aSHong Zhang /*@
1743b5b4017aSHong Zhang   TSForwardSetInitialSensitivities - Set initial values for tangent linear sensitivities
1744b5b4017aSHong Zhang 
1745c3339decSBarry Smith   Collective
1746b5b4017aSHong Zhang 
1747d8d19677SJose E. Roman   Input Parameters:
1748bcf0153eSBarry Smith + ts - the `TS` context obtained from `TSCreate()`
1749b5b4017aSHong Zhang - didp - parametric sensitivities of the initial condition
1750b5b4017aSHong Zhang 
1751b5b4017aSHong Zhang   Level: intermediate
1752b5b4017aSHong Zhang 
1753bcf0153eSBarry Smith   Notes:
1754bcf0153eSBarry Smith   `TSSolve()` allows users to pass the initial solution directly to `TS`. But the tangent linear variables cannot be initialized in this way.
1755bcf0153eSBarry Smith    This function is used to set initial values for tangent linear variables.
1756b5b4017aSHong Zhang 
1757bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSForwardSetSensitivities()`
1758b5b4017aSHong Zhang @*/
1759d71ae5a4SJacob Faibussowitsch PetscErrorCode TSForwardSetInitialSensitivities(TS ts, Mat didp)
1760d71ae5a4SJacob Faibussowitsch {
1761362febeeSStefano Zampini   PetscFunctionBegin;
1762b5b4017aSHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1763b5b4017aSHong Zhang   PetscValidHeaderSpecific(didp, MAT_CLASSID, 2);
176448a46eb9SPierre Jolivet   if (!ts->mat_sensip) PetscCall(TSForwardSetSensitivities(ts, PETSC_DEFAULT, didp));
17653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1766b5b4017aSHong Zhang }
17676affb6f8SHong Zhang 
17686affb6f8SHong Zhang /*@
17696affb6f8SHong Zhang    TSForwardGetStages - Get the number of stages and the tangent linear sensitivities at the intermediate stages
17706affb6f8SHong Zhang 
17716affb6f8SHong Zhang    Input Parameter:
1772bcf0153eSBarry Smith .  ts - the `TS` context obtained from `TSCreate()`
17736affb6f8SHong Zhang 
17746affb6f8SHong Zhang    Output Parameters:
1775cd4cee2dSHong Zhang +  ns - number of stages
17766affb6f8SHong Zhang -  S - tangent linear sensitivities at the intermediate stages
17776affb6f8SHong Zhang 
17786affb6f8SHong Zhang    Level: advanced
17796affb6f8SHong Zhang 
1780bcf0153eSBarry Smith .seealso: `TS`
17816affb6f8SHong Zhang @*/
1782d71ae5a4SJacob Faibussowitsch PetscErrorCode TSForwardGetStages(TS ts, PetscInt *ns, Mat **S)
1783d71ae5a4SJacob Faibussowitsch {
17846affb6f8SHong Zhang   PetscFunctionBegin;
17856affb6f8SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
17866affb6f8SHong Zhang 
17876affb6f8SHong Zhang   if (!ts->ops->getstages) *S = NULL;
1788dbbe0bcdSBarry Smith   else PetscUseTypeMethod(ts, forwardgetstages, ns, S);
17893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
17906affb6f8SHong Zhang }
1791cd4cee2dSHong Zhang 
1792cd4cee2dSHong Zhang /*@
1793bcf0153eSBarry Smith    TSCreateQuadratureTS - Create a sub-`TS` that evaluates integrals over time
1794cd4cee2dSHong Zhang 
1795d8d19677SJose E. Roman    Input Parameters:
1796bcf0153eSBarry Smith +  ts - the `TS` context obtained from `TSCreate()`
1797cd4cee2dSHong Zhang -  fwd - flag indicating whether to evaluate cost integral in the forward run or the adjoint run
1798cd4cee2dSHong Zhang 
1799cd4cee2dSHong Zhang    Output Parameters:
1800bcf0153eSBarry Smith .  quadts - the child `TS` context
1801cd4cee2dSHong Zhang 
1802cd4cee2dSHong Zhang    Level: intermediate
1803cd4cee2dSHong Zhang 
1804bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSGetQuadratureTS()`
1805cd4cee2dSHong Zhang @*/
1806d71ae5a4SJacob Faibussowitsch PetscErrorCode TSCreateQuadratureTS(TS ts, PetscBool fwd, TS *quadts)
1807d71ae5a4SJacob Faibussowitsch {
1808cd4cee2dSHong Zhang   char prefix[128];
1809cd4cee2dSHong Zhang 
1810cd4cee2dSHong Zhang   PetscFunctionBegin;
1811cd4cee2dSHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1812064a246eSJacob Faibussowitsch   PetscValidPointer(quadts, 3);
18139566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts->quadraturets));
18149566063dSJacob Faibussowitsch   PetscCall(TSCreate(PetscObjectComm((PetscObject)ts), &ts->quadraturets));
18159566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)ts->quadraturets, (PetscObject)ts, 1));
18169566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(prefix, sizeof(prefix), "%squad_", ((PetscObject)ts)->prefix ? ((PetscObject)ts)->prefix : ""));
18179566063dSJacob Faibussowitsch   PetscCall(TSSetOptionsPrefix(ts->quadraturets, prefix));
1818cd4cee2dSHong Zhang   *quadts = ts->quadraturets;
1819cd4cee2dSHong Zhang 
1820cd4cee2dSHong Zhang   if (ts->numcost) {
18219566063dSJacob Faibussowitsch     PetscCall(VecCreateSeq(PETSC_COMM_SELF, ts->numcost, &(*quadts)->vec_sol));
1822cd4cee2dSHong Zhang   } else {
18239566063dSJacob Faibussowitsch     PetscCall(VecCreateSeq(PETSC_COMM_SELF, 1, &(*quadts)->vec_sol));
1824cd4cee2dSHong Zhang   }
1825cd4cee2dSHong Zhang   ts->costintegralfwd = fwd;
18263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1827cd4cee2dSHong Zhang }
1828cd4cee2dSHong Zhang 
1829cd4cee2dSHong Zhang /*@
1830bcf0153eSBarry Smith    TSGetQuadratureTS - Return the sub-`TS` that evaluates integrals over time
1831cd4cee2dSHong Zhang 
1832cd4cee2dSHong Zhang    Input Parameter:
1833bcf0153eSBarry Smith .  ts - the `TS` context obtained from `TSCreate()`
1834cd4cee2dSHong Zhang 
1835cd4cee2dSHong Zhang    Output Parameters:
1836cd4cee2dSHong Zhang +  fwd - flag indicating whether to evaluate cost integral in the forward run or the adjoint run
1837bcf0153eSBarry Smith -  quadts - the child `TS` context
1838cd4cee2dSHong Zhang 
1839cd4cee2dSHong Zhang    Level: intermediate
1840cd4cee2dSHong Zhang 
1841bcf0153eSBarry Smith .seealso: [](chapter_ts), `TSCreateQuadratureTS()`
1842cd4cee2dSHong Zhang @*/
1843d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGetQuadratureTS(TS ts, PetscBool *fwd, TS *quadts)
1844d71ae5a4SJacob Faibussowitsch {
1845cd4cee2dSHong Zhang   PetscFunctionBegin;
1846cd4cee2dSHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1847cd4cee2dSHong Zhang   if (fwd) *fwd = ts->costintegralfwd;
1848cd4cee2dSHong Zhang   if (quadts) *quadts = ts->quadraturets;
18493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1850cd4cee2dSHong Zhang }
1851ba0675f6SHong Zhang 
1852ba0675f6SHong Zhang /*@
1853bcf0153eSBarry Smith    TSComputeSNESJacobian - Compute the Jacobian needed for the `SNESSolve()` in `TS`
1854bcf0153eSBarry Smith 
1855c3339decSBarry Smith    Collective
1856ba0675f6SHong Zhang 
1857ba0675f6SHong Zhang    Input Parameters:
1858bcf0153eSBarry Smith +  ts - the `TS` context obtained from `TSCreate()`
1859ba0675f6SHong Zhang -  x - state vector
1860ba0675f6SHong Zhang 
1861ba0675f6SHong Zhang    Output Parameters:
1862ba0675f6SHong Zhang +  J - Jacobian matrix
1863ba0675f6SHong Zhang -  Jpre - preconditioning matrix for J (may be same as J)
1864ba0675f6SHong Zhang 
1865ba0675f6SHong Zhang    Level: developer
1866ba0675f6SHong Zhang 
1867bcf0153eSBarry Smith    Note:
1868bcf0153eSBarry Smith    Uses finite differencing when `TS` Jacobian is not available.
1869bcf0153eSBarry Smith 
1870bcf0153eSBarry Smith .seealso: `SNES`, `TS`, `SNESSetJacobian()`, TSSetRHSJacobian()`, `TSSetIJacobian()`
1871ba0675f6SHong Zhang @*/
1872d71ae5a4SJacob Faibussowitsch PetscErrorCode TSComputeSNESJacobian(TS ts, Vec x, Mat J, Mat Jpre)
1873d71ae5a4SJacob Faibussowitsch {
1874ba0675f6SHong Zhang   SNES snes                                          = ts->snes;
1875ba0675f6SHong Zhang   PetscErrorCode (*jac)(SNES, Vec, Mat, Mat, void *) = NULL;
1876ba0675f6SHong Zhang 
1877ba0675f6SHong Zhang   PetscFunctionBegin;
1878ba0675f6SHong Zhang   /*
1879ba0675f6SHong Zhang     Unlike implicit methods, explicit methods do not have SNESMatFDColoring in the snes object
1880ba0675f6SHong Zhang     because SNESSolve() has not been called yet; so querying SNESMatFDColoring does not work for
1881da81f932SPierre Jolivet     explicit methods. Instead, we check the Jacobian compute function directly to determine if FD
1882ba0675f6SHong Zhang     coloring is used.
1883ba0675f6SHong Zhang   */
18849566063dSJacob Faibussowitsch   PetscCall(SNESGetJacobian(snes, NULL, NULL, &jac, NULL));
1885ba0675f6SHong Zhang   if (jac == SNESComputeJacobianDefaultColor) {
1886ba0675f6SHong Zhang     Vec f;
18879566063dSJacob Faibussowitsch     PetscCall(SNESSetSolution(snes, x));
18889566063dSJacob Faibussowitsch     PetscCall(SNESGetFunction(snes, &f, NULL, NULL));
1889ba0675f6SHong Zhang     /* Force MatFDColoringApply to evaluate the SNES residual function for the base vector */
18909566063dSJacob Faibussowitsch     PetscCall(SNESComputeFunction(snes, x, f));
1891ba0675f6SHong Zhang   }
18929566063dSJacob Faibussowitsch   PetscCall(SNESComputeJacobian(snes, x, J, Jpre));
18933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1894ba0675f6SHong Zhang }
1895