xref: /petsc/src/ts/interface/sensitivity/tssen.c (revision ba38deedd2838b0f18e3744cb0bafd6392690fde)
114d0ab18SJacob Faibussowitsch #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   Level: intermediate
22814e85d6SHong Zhang 
23bcf0153eSBarry Smith   Note:
2414d0ab18SJacob Faibussowitsch   `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`
25814e85d6SHong Zhang 
2614d0ab18SJacob Faibussowitsch .seealso: [](ch_ts), `TS`, `TSRHSJacobianP`, `TSGetRHSJacobianP()`
27814e85d6SHong Zhang @*/
2814d0ab18SJacob Faibussowitsch PetscErrorCode TSSetRHSJacobianP(TS ts, Mat Amat, TSRHSJacobianP func, void *ctx)
29d71ae5a4SJacob Faibussowitsch {
30814e85d6SHong Zhang   PetscFunctionBegin;
31814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
32814e85d6SHong Zhang   PetscValidHeaderSpecific(Amat, MAT_CLASSID, 2);
33814e85d6SHong Zhang 
34814e85d6SHong Zhang   ts->rhsjacobianp    = func;
35814e85d6SHong Zhang   ts->rhsjacobianpctx = ctx;
36814e85d6SHong Zhang   if (Amat) {
379566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)Amat));
389566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&ts->Jacprhs));
3933f72c5dSHong Zhang     ts->Jacprhs = Amat;
40814e85d6SHong Zhang   }
413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
42814e85d6SHong Zhang }
43814e85d6SHong Zhang 
44814e85d6SHong Zhang /*@C
45cd4cee2dSHong 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.
46cd4cee2dSHong Zhang 
47c3339decSBarry Smith   Logically Collective
48cd4cee2dSHong Zhang 
49f899ff85SJose E. Roman   Input Parameter:
50bcf0153eSBarry Smith . ts - `TS` context obtained from `TSCreate()`
51cd4cee2dSHong Zhang 
52cd4cee2dSHong Zhang   Output Parameters:
53cd4cee2dSHong Zhang + Amat - JacobianP matrix
54cd4cee2dSHong Zhang . func - function
55cd4cee2dSHong Zhang - ctx  - [optional] user-defined function context
56cd4cee2dSHong Zhang 
57cd4cee2dSHong Zhang   Level: intermediate
58cd4cee2dSHong Zhang 
59bcf0153eSBarry Smith   Note:
6014d0ab18SJacob Faibussowitsch   `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`
61cd4cee2dSHong Zhang 
6242747ad1SJacob Faibussowitsch .seealso: [](ch_ts), `TSSetRHSJacobianP()`, `TS`, `TSRHSJacobianP`
63cd4cee2dSHong Zhang @*/
6414d0ab18SJacob Faibussowitsch PetscErrorCode TSGetRHSJacobianP(TS ts, Mat *Amat, TSRHSJacobianP *func, void **ctx)
65d71ae5a4SJacob Faibussowitsch {
66cd4cee2dSHong Zhang   PetscFunctionBegin;
67cd4cee2dSHong Zhang   if (func) *func = ts->rhsjacobianp;
68cd4cee2dSHong Zhang   if (ctx) *ctx = ts->rhsjacobianpctx;
69cd4cee2dSHong Zhang   if (Amat) *Amat = ts->Jacprhs;
703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
71cd4cee2dSHong Zhang }
72cd4cee2dSHong Zhang 
73cd4cee2dSHong Zhang /*@C
74814e85d6SHong Zhang   TSComputeRHSJacobianP - Runs the user-defined JacobianP function.
75814e85d6SHong Zhang 
76c3339decSBarry Smith   Collective
77814e85d6SHong Zhang 
78814e85d6SHong Zhang   Input Parameters:
792fe279fdSBarry Smith + ts - The `TS` context obtained from `TSCreate()`
802fe279fdSBarry Smith . t  - the time
812fe279fdSBarry Smith - U  - the solution at which to compute the Jacobian
822fe279fdSBarry Smith 
832fe279fdSBarry Smith   Output Parameter:
842fe279fdSBarry Smith . Amat - the computed Jacobian
85814e85d6SHong Zhang 
86814e85d6SHong Zhang   Level: developer
87814e85d6SHong Zhang 
881cc06b55SBarry Smith .seealso: [](ch_ts), `TSSetRHSJacobianP()`, `TS`
89814e85d6SHong Zhang @*/
90d71ae5a4SJacob Faibussowitsch PetscErrorCode TSComputeRHSJacobianP(TS ts, PetscReal t, Vec U, Mat Amat)
91d71ae5a4SJacob Faibussowitsch {
92814e85d6SHong Zhang   PetscFunctionBegin;
933ba16761SJacob Faibussowitsch   if (!Amat) PetscFunctionReturn(PETSC_SUCCESS);
94814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
95c9ad9de0SHong Zhang   PetscValidHeaderSpecific(U, VEC_CLASSID, 3);
96814e85d6SHong Zhang 
9780ab5e31SHong Zhang   if (ts->rhsjacobianp) PetscCallBack("TS callback JacobianP for sensitivity analysis", (*ts->rhsjacobianp)(ts, t, U, Amat, ts->rhsjacobianpctx));
9880ab5e31SHong Zhang   else {
9980ab5e31SHong Zhang     PetscBool assembled;
10080ab5e31SHong Zhang     PetscCall(MatZeroEntries(Amat));
10180ab5e31SHong Zhang     PetscCall(MatAssembled(Amat, &assembled));
10280ab5e31SHong Zhang     if (!assembled) {
10380ab5e31SHong Zhang       PetscCall(MatAssemblyBegin(Amat, MAT_FINAL_ASSEMBLY));
10480ab5e31SHong Zhang       PetscCall(MatAssemblyEnd(Amat, MAT_FINAL_ASSEMBLY));
10580ab5e31SHong Zhang     }
10680ab5e31SHong Zhang   }
1073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
108814e85d6SHong Zhang }
109814e85d6SHong Zhang 
110814e85d6SHong Zhang /*@C
11133f72c5dSHong 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.
11233f72c5dSHong Zhang 
113c3339decSBarry Smith   Logically Collective
11433f72c5dSHong Zhang 
11533f72c5dSHong Zhang   Input Parameters:
116bcf0153eSBarry Smith + ts   - `TS` context obtained from `TSCreate()`
11733f72c5dSHong Zhang . Amat - JacobianP matrix
11833f72c5dSHong Zhang . func - function
11933f72c5dSHong Zhang - ctx  - [optional] user-defined function context
12033f72c5dSHong Zhang 
12120f4b53cSBarry Smith   Calling sequence of `func`:
12214d0ab18SJacob Faibussowitsch + ts    - the `TS` context
12314d0ab18SJacob Faibussowitsch . t     - current timestep
12433f72c5dSHong Zhang . U     - input vector (current ODE solution)
12533f72c5dSHong Zhang . Udot  - time derivative of state vector
12633f72c5dSHong Zhang . shift - shift to apply, see note below
12733f72c5dSHong Zhang . A     - output matrix
12833f72c5dSHong Zhang - ctx   - [optional] user-defined function context
12933f72c5dSHong Zhang 
13033f72c5dSHong Zhang   Level: intermediate
13133f72c5dSHong Zhang 
132bcf0153eSBarry Smith   Note:
13333f72c5dSHong 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
13433f72c5dSHong Zhang 
1351cc06b55SBarry Smith .seealso: [](ch_ts), `TSSetRHSJacobianP()`, `TS`
13633f72c5dSHong Zhang @*/
13714d0ab18SJacob Faibussowitsch PetscErrorCode TSSetIJacobianP(TS ts, Mat Amat, PetscErrorCode (*func)(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal shift, Mat A, void *ctx), void *ctx)
138d71ae5a4SJacob Faibussowitsch {
13933f72c5dSHong Zhang   PetscFunctionBegin;
14033f72c5dSHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
14133f72c5dSHong Zhang   PetscValidHeaderSpecific(Amat, MAT_CLASSID, 2);
14233f72c5dSHong Zhang 
14333f72c5dSHong Zhang   ts->ijacobianp    = func;
14433f72c5dSHong Zhang   ts->ijacobianpctx = ctx;
14533f72c5dSHong Zhang   if (Amat) {
1469566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)Amat));
1479566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&ts->Jacp));
14833f72c5dSHong Zhang     ts->Jacp = Amat;
14933f72c5dSHong Zhang   }
1503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15133f72c5dSHong Zhang }
15233f72c5dSHong Zhang 
15333f72c5dSHong Zhang /*@C
15433f72c5dSHong Zhang   TSComputeIJacobianP - Runs the user-defined IJacobianP function.
15533f72c5dSHong Zhang 
156c3339decSBarry Smith   Collective
15733f72c5dSHong Zhang 
15833f72c5dSHong Zhang   Input Parameters:
159bcf0153eSBarry Smith + ts    - the `TS` context
16033f72c5dSHong Zhang . t     - current timestep
16133f72c5dSHong Zhang . U     - state vector
16233f72c5dSHong Zhang . Udot  - time derivative of state vector
16333f72c5dSHong Zhang . shift - shift to apply, see note below
16480ab5e31SHong Zhang - imex  - flag indicates if the method is IMEX so that the RHSJacobianP should be kept separate
16533f72c5dSHong Zhang 
1662fe279fdSBarry Smith   Output Parameter:
167b43aa488SJacob Faibussowitsch . Amat - Jacobian matrix
16833f72c5dSHong Zhang 
16933f72c5dSHong Zhang   Level: developer
17033f72c5dSHong Zhang 
1711cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSetIJacobianP()`
17233f72c5dSHong Zhang @*/
173d71ae5a4SJacob Faibussowitsch PetscErrorCode TSComputeIJacobianP(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal shift, Mat Amat, PetscBool imex)
174d71ae5a4SJacob Faibussowitsch {
17533f72c5dSHong Zhang   PetscFunctionBegin;
1763ba16761SJacob Faibussowitsch   if (!Amat) PetscFunctionReturn(PETSC_SUCCESS);
17733f72c5dSHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
17833f72c5dSHong Zhang   PetscValidHeaderSpecific(U, VEC_CLASSID, 3);
17933f72c5dSHong Zhang   PetscValidHeaderSpecific(Udot, VEC_CLASSID, 4);
18033f72c5dSHong Zhang 
1819566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(TS_JacobianPEval, ts, U, Amat, 0));
18248a46eb9SPierre Jolivet   if (ts->ijacobianp) PetscCallBack("TS callback JacobianP for sensitivity analysis", (*ts->ijacobianp)(ts, t, U, Udot, shift, Amat, ts->ijacobianpctx));
18333f72c5dSHong Zhang   if (imex) {
18433f72c5dSHong Zhang     if (!ts->ijacobianp) { /* system was written as Udot = G(t,U) */
18533f72c5dSHong Zhang       PetscBool assembled;
1869566063dSJacob Faibussowitsch       PetscCall(MatZeroEntries(Amat));
1879566063dSJacob Faibussowitsch       PetscCall(MatAssembled(Amat, &assembled));
18833f72c5dSHong Zhang       if (!assembled) {
1899566063dSJacob Faibussowitsch         PetscCall(MatAssemblyBegin(Amat, MAT_FINAL_ASSEMBLY));
1909566063dSJacob Faibussowitsch         PetscCall(MatAssemblyEnd(Amat, MAT_FINAL_ASSEMBLY));
19133f72c5dSHong Zhang       }
19233f72c5dSHong Zhang     }
19333f72c5dSHong Zhang   } else {
1941baa6e33SBarry Smith     if (ts->rhsjacobianp) PetscCall(TSComputeRHSJacobianP(ts, t, U, ts->Jacprhs));
19533f72c5dSHong Zhang     if (ts->Jacprhs == Amat) { /* No IJacobian, so we only have the RHS matrix */
1969566063dSJacob Faibussowitsch       PetscCall(MatScale(Amat, -1));
19733f72c5dSHong Zhang     } else if (ts->Jacprhs) { /* Both IJacobian and RHSJacobian */
19833f72c5dSHong Zhang       MatStructure axpy = DIFFERENT_NONZERO_PATTERN;
19933f72c5dSHong Zhang       if (!ts->ijacobianp) { /* No IJacobianp provided, but we have a separate RHS matrix */
2009566063dSJacob Faibussowitsch         PetscCall(MatZeroEntries(Amat));
20133f72c5dSHong Zhang       }
2029566063dSJacob Faibussowitsch       PetscCall(MatAXPY(Amat, -1, ts->Jacprhs, axpy));
20333f72c5dSHong Zhang     }
20433f72c5dSHong Zhang   }
2059566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(TS_JacobianPEval, ts, U, Amat, 0));
2063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
20733f72c5dSHong Zhang }
20833f72c5dSHong Zhang 
20933f72c5dSHong Zhang /*@C
210814e85d6SHong Zhang   TSSetCostIntegrand - Sets the routine for evaluating the integral term in one or more cost functions
211814e85d6SHong Zhang 
212c3339decSBarry Smith   Logically Collective
213814e85d6SHong Zhang 
214814e85d6SHong Zhang   Input Parameters:
215bcf0153eSBarry Smith + ts           - the `TS` context obtained from `TSCreate()`
216814e85d6SHong Zhang . numcost      - number of gradients to be computed, this is the number of cost functions
217814e85d6SHong Zhang . costintegral - vector that stores the integral values
218814e85d6SHong Zhang . rf           - routine for evaluating the integrand function
219c9ad9de0SHong Zhang . drduf        - function that computes the gradients of the r's with respect to u
2202fe279fdSBarry Smith . 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`)
221814e85d6SHong Zhang . fwd          - flag indicating whether to evaluate cost integral in the forward run or the adjoint run
2222fe279fdSBarry Smith - ctx          - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`)
223814e85d6SHong Zhang 
22420f4b53cSBarry Smith   Calling sequence of `rf`:
22520f4b53cSBarry Smith $   PetscErrorCode rf(TS ts, PetscReal t, Vec U, Vec F, oid *ctx)
226814e85d6SHong Zhang 
22720f4b53cSBarry Smith   Calling sequence of `drduf`:
22820f4b53cSBarry Smith $   PetscErroCode drduf(TS ts, PetscReal t, Vec U, Vec *dRdU, void *ctx)
229814e85d6SHong Zhang 
23020f4b53cSBarry Smith   Calling sequence of `drdpf`:
23120f4b53cSBarry Smith $   PetscErroCode drdpf(TS ts, PetscReal t, Vec U, Vec *dRdP, void *ctx)
232814e85d6SHong Zhang 
233cd4cee2dSHong Zhang   Level: deprecated
234814e85d6SHong Zhang 
235bcf0153eSBarry Smith   Note:
236814e85d6SHong Zhang   For optimization there is usually a single cost function (numcost = 1). For sensitivities there may be multiple cost functions
237814e85d6SHong Zhang 
2381cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSetRHSJacobianP()`, `TSGetCostGradients()`, `TSSetCostGradients()`
239814e85d6SHong Zhang @*/
240d71ae5a4SJacob 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)
241d71ae5a4SJacob Faibussowitsch {
242814e85d6SHong Zhang   PetscFunctionBegin;
243814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
244814e85d6SHong Zhang   if (costintegral) PetscValidHeaderSpecific(costintegral, VEC_CLASSID, 3);
2453c633725SBarry 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()");
246814e85d6SHong Zhang   if (!ts->numcost) ts->numcost = numcost;
247814e85d6SHong Zhang 
248814e85d6SHong Zhang   if (costintegral) {
2499566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)costintegral));
2509566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&ts->vec_costintegral));
251814e85d6SHong Zhang     ts->vec_costintegral = costintegral;
252814e85d6SHong Zhang   } else {
253814e85d6SHong Zhang     if (!ts->vec_costintegral) { /* Create a seq vec if user does not provide one */
2549566063dSJacob Faibussowitsch       PetscCall(VecCreateSeq(PETSC_COMM_SELF, numcost, &ts->vec_costintegral));
255814e85d6SHong Zhang     } else {
2569566063dSJacob Faibussowitsch       PetscCall(VecSet(ts->vec_costintegral, 0.0));
257814e85d6SHong Zhang     }
258814e85d6SHong Zhang   }
259814e85d6SHong Zhang   if (!ts->vec_costintegrand) {
2609566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(ts->vec_costintegral, &ts->vec_costintegrand));
261814e85d6SHong Zhang   } else {
2629566063dSJacob Faibussowitsch     PetscCall(VecSet(ts->vec_costintegrand, 0.0));
263814e85d6SHong Zhang   }
264814e85d6SHong Zhang   ts->costintegralfwd  = fwd; /* Evaluate the cost integral in forward run if fwd is true */
265814e85d6SHong Zhang   ts->costintegrand    = rf;
266814e85d6SHong Zhang   ts->costintegrandctx = ctx;
267c9ad9de0SHong Zhang   ts->drdufunction     = drduf;
268814e85d6SHong Zhang   ts->drdpfunction     = drdpf;
2693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
270814e85d6SHong Zhang }
271814e85d6SHong Zhang 
272b5b4017aSHong Zhang /*@C
273814e85d6SHong Zhang   TSGetCostIntegral - Returns the values of the integral term in the cost functions.
274814e85d6SHong Zhang   It is valid to call the routine after a backward run.
275814e85d6SHong Zhang 
276814e85d6SHong Zhang   Not Collective
277814e85d6SHong Zhang 
278814e85d6SHong Zhang   Input Parameter:
279bcf0153eSBarry Smith . ts - the `TS` context obtained from `TSCreate()`
280814e85d6SHong Zhang 
281814e85d6SHong Zhang   Output Parameter:
282814e85d6SHong Zhang . v - the vector containing the integrals for each cost function
283814e85d6SHong Zhang 
284814e85d6SHong Zhang   Level: intermediate
285814e85d6SHong Zhang 
2861cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSAdjointSolve()`, ``TSSetCostIntegrand()`
287814e85d6SHong Zhang @*/
288d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGetCostIntegral(TS ts, Vec *v)
289d71ae5a4SJacob Faibussowitsch {
290cd4cee2dSHong Zhang   TS quadts;
291cd4cee2dSHong Zhang 
292814e85d6SHong Zhang   PetscFunctionBegin;
293814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
2944f572ea9SToby Isaac   PetscAssertPointer(v, 2);
2959566063dSJacob Faibussowitsch   PetscCall(TSGetQuadratureTS(ts, NULL, &quadts));
296cd4cee2dSHong Zhang   *v = quadts->vec_sol;
2973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
298814e85d6SHong Zhang }
299814e85d6SHong Zhang 
300b5b4017aSHong Zhang /*@C
301814e85d6SHong Zhang   TSComputeCostIntegrand - Evaluates the integral function in the cost functions.
302814e85d6SHong Zhang 
303814e85d6SHong Zhang   Input Parameters:
304bcf0153eSBarry Smith + ts - the `TS` context
305814e85d6SHong Zhang . t  - current time
306c9ad9de0SHong Zhang - U  - state vector, i.e. current solution
307814e85d6SHong Zhang 
308814e85d6SHong Zhang   Output Parameter:
309c9ad9de0SHong Zhang . Q - vector of size numcost to hold the outputs
310814e85d6SHong Zhang 
311bcf0153eSBarry Smith   Level: deprecated
312bcf0153eSBarry Smith 
313bcf0153eSBarry Smith   Note:
314814e85d6SHong Zhang   Most users should not need to explicitly call this routine, as it
315814e85d6SHong Zhang   is used internally within the sensitivity analysis context.
316814e85d6SHong Zhang 
3171cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSAdjointSolve()`, `TSSetCostIntegrand()`
318814e85d6SHong Zhang @*/
319d71ae5a4SJacob Faibussowitsch PetscErrorCode TSComputeCostIntegrand(TS ts, PetscReal t, Vec U, Vec Q)
320d71ae5a4SJacob Faibussowitsch {
321814e85d6SHong Zhang   PetscFunctionBegin;
322814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
323c9ad9de0SHong Zhang   PetscValidHeaderSpecific(U, VEC_CLASSID, 3);
324c9ad9de0SHong Zhang   PetscValidHeaderSpecific(Q, VEC_CLASSID, 4);
325814e85d6SHong Zhang 
3269566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(TS_FunctionEval, ts, U, Q, 0));
327792fecdfSBarry Smith   if (ts->costintegrand) PetscCallBack("TS callback integrand in the cost function", (*ts->costintegrand)(ts, t, U, Q, ts->costintegrandctx));
328ef1023bdSBarry Smith   else PetscCall(VecZeroEntries(Q));
3299566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(TS_FunctionEval, ts, U, Q, 0));
3303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
331814e85d6SHong Zhang }
332814e85d6SHong Zhang 
33314d0ab18SJacob Faibussowitsch // PetscClangLinter pragma disable: -fdoc-*
334b5b4017aSHong Zhang /*@C
335bcf0153eSBarry Smith   TSComputeDRDUFunction - Deprecated, use `TSGetQuadratureTS()` then `TSComputeRHSJacobian()`
336814e85d6SHong Zhang 
337d76be551SHong Zhang   Level: deprecated
338814e85d6SHong Zhang 
339814e85d6SHong Zhang @*/
340d71ae5a4SJacob Faibussowitsch PetscErrorCode TSComputeDRDUFunction(TS ts, PetscReal t, Vec U, Vec *DRDU)
341d71ae5a4SJacob Faibussowitsch {
342814e85d6SHong Zhang   PetscFunctionBegin;
3433ba16761SJacob Faibussowitsch   if (!DRDU) PetscFunctionReturn(PETSC_SUCCESS);
344814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
345c9ad9de0SHong Zhang   PetscValidHeaderSpecific(U, VEC_CLASSID, 3);
346814e85d6SHong Zhang 
347792fecdfSBarry Smith   PetscCallBack("TS callback DRDU for sensitivity analysis", (*ts->drdufunction)(ts, t, U, DRDU, ts->costintegrandctx));
3483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
349814e85d6SHong Zhang }
350814e85d6SHong Zhang 
35114d0ab18SJacob Faibussowitsch // PetscClangLinter pragma disable: -fdoc-*
352b5b4017aSHong Zhang /*@C
353bcf0153eSBarry Smith   TSComputeDRDPFunction - Deprecated, use `TSGetQuadratureTS()` then `TSComputeRHSJacobianP()`
354814e85d6SHong Zhang 
355d76be551SHong Zhang   Level: deprecated
356814e85d6SHong Zhang 
357814e85d6SHong Zhang @*/
358d71ae5a4SJacob Faibussowitsch PetscErrorCode TSComputeDRDPFunction(TS ts, PetscReal t, Vec U, Vec *DRDP)
359d71ae5a4SJacob Faibussowitsch {
360814e85d6SHong Zhang   PetscFunctionBegin;
3613ba16761SJacob Faibussowitsch   if (!DRDP) PetscFunctionReturn(PETSC_SUCCESS);
362814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
363c9ad9de0SHong Zhang   PetscValidHeaderSpecific(U, VEC_CLASSID, 3);
364814e85d6SHong Zhang 
365792fecdfSBarry Smith   PetscCallBack("TS callback DRDP for sensitivity analysis", (*ts->drdpfunction)(ts, t, U, DRDP, ts->costintegrandctx));
3663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
367814e85d6SHong Zhang }
368814e85d6SHong Zhang 
36914d0ab18SJacob Faibussowitsch // PetscClangLinter pragma disable: -fdoc-param-list-func-parameter-documentation
37014d0ab18SJacob Faibussowitsch // PetscClangLinter pragma disable: -fdoc-section-header-unknown
371b5b4017aSHong Zhang /*@C
372b5b4017aSHong 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.
373b5b4017aSHong Zhang 
374c3339decSBarry Smith   Logically Collective
375b5b4017aSHong Zhang 
376b5b4017aSHong Zhang   Input Parameters:
377bcf0153eSBarry Smith + ts   - `TS` context obtained from `TSCreate()`
378b5b4017aSHong Zhang . ihp1 - an array of vectors storing the result of vector-Hessian-vector product for F_UU
379b5b4017aSHong Zhang . hessianproductfunc1 - vector-Hessian-vector product function for F_UU
380b5b4017aSHong Zhang . ihp2 - an array of vectors storing the result of vector-Hessian-vector product for F_UP
381b5b4017aSHong Zhang . hessianproductfunc2 - vector-Hessian-vector product function for F_UP
382b5b4017aSHong Zhang . ihp3 - an array of vectors storing the result of vector-Hessian-vector product for F_PU
383b5b4017aSHong Zhang . hessianproductfunc3 - vector-Hessian-vector product function for F_PU
384b5b4017aSHong Zhang . ihp4 - an array of vectors storing the result of vector-Hessian-vector product for F_PP
385f0fc11ceSJed Brown - hessianproductfunc4 - vector-Hessian-vector product function for F_PP
386b5b4017aSHong Zhang 
38714d0ab18SJacob Faibussowitsch   Calling sequence of `ihessianproductfunc1`:
38814d0ab18SJacob Faibussowitsch + ts  - the `TS` context
389b5b4017aSHong Zhang + t   - current timestep
390b5b4017aSHong Zhang . U   - input vector (current ODE solution)
391369cf35fSHong Zhang . Vl  - an array of input vectors to be left-multiplied with the Hessian
392b5b4017aSHong Zhang . Vr  - input vector to be right-multiplied with the Hessian
393369cf35fSHong Zhang . VHV - an array of output vectors for vector-Hessian-vector product
394b5b4017aSHong Zhang - ctx - [optional] user-defined function context
395b5b4017aSHong Zhang 
396b5b4017aSHong Zhang   Level: intermediate
397b5b4017aSHong Zhang 
398369cf35fSHong Zhang   Notes:
39914d0ab18SJacob Faibussowitsch   All other functions have the same calling sequence as `rhhessianproductfunc1`, so their
40014d0ab18SJacob Faibussowitsch   descriptions are omitted for brevity.
40114d0ab18SJacob Faibussowitsch 
402369cf35fSHong Zhang   The first Hessian function and the working array are required.
403369cf35fSHong Zhang   As an example to implement the callback functions, the second callback function calculates the vector-Hessian-vector product
404369cf35fSHong Zhang   $ Vl_n^T*F_UP*Vr
405369cf35fSHong 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.
406369cf35fSHong Zhang   Each entry of F_UP corresponds to the derivative
407369cf35fSHong Zhang   $ F_UP[i][j][k] = \frac{\partial^2 F[i]}{\partial U[j] \partial P[k]}.
408369cf35fSHong 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
409369cf35fSHong Zhang   $ VHV_n[j] = \sum_i \sum_k {Vl_n[i] * F_UP[i][j][k] * Vr[k]}
410369cf35fSHong Zhang   If the cost function is a scalar, there will be only one vector in Vl and VHV.
411b5b4017aSHong Zhang 
4121cc06b55SBarry Smith .seealso: [](ch_ts), `TS`
413b5b4017aSHong Zhang @*/
41414d0ab18SJacob Faibussowitsch PetscErrorCode TSSetIHessianProduct(TS ts, Vec *ihp1, PetscErrorCode (*ihessianproductfunc1)(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV, void *ctx), 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)
415d71ae5a4SJacob Faibussowitsch {
416b5b4017aSHong Zhang   PetscFunctionBegin;
417b5b4017aSHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
4184f572ea9SToby Isaac   PetscAssertPointer(ihp1, 2);
419b5b4017aSHong Zhang 
420b5b4017aSHong Zhang   ts->ihessianproductctx = ctx;
421b5b4017aSHong Zhang   if (ihp1) ts->vecs_fuu = ihp1;
422b5b4017aSHong Zhang   if (ihp2) ts->vecs_fup = ihp2;
423b5b4017aSHong Zhang   if (ihp3) ts->vecs_fpu = ihp3;
424b5b4017aSHong Zhang   if (ihp4) ts->vecs_fpp = ihp4;
425b5b4017aSHong Zhang   ts->ihessianproduct_fuu = ihessianproductfunc1;
426b5b4017aSHong Zhang   ts->ihessianproduct_fup = ihessianproductfunc2;
427b5b4017aSHong Zhang   ts->ihessianproduct_fpu = ihessianproductfunc3;
428b5b4017aSHong Zhang   ts->ihessianproduct_fpp = ihessianproductfunc4;
4293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
430b5b4017aSHong Zhang }
431b5b4017aSHong Zhang 
432b5b4017aSHong Zhang /*@C
433b1e111ebSHong Zhang   TSComputeIHessianProductFunctionUU - Runs the user-defined vector-Hessian-vector product function for Fuu.
434b5b4017aSHong Zhang 
435c3339decSBarry Smith   Collective
436b5b4017aSHong Zhang 
437b5b4017aSHong Zhang   Input Parameters:
4382fe279fdSBarry Smith + ts - The `TS` context obtained from `TSCreate()`
4392fe279fdSBarry Smith . t  - the time
4402fe279fdSBarry Smith . U  - the solution at which to compute the Hessian product
4412fe279fdSBarry Smith . Vl - the array of input vectors to be multiplied with the Hessian from the left
4422fe279fdSBarry Smith - Vr - the input vector to be multiplied with the Hessian from the right
4432fe279fdSBarry Smith 
4442fe279fdSBarry Smith   Output Parameter:
445b43aa488SJacob Faibussowitsch . VHV - the array of output vectors that store the Hessian product
446b5b4017aSHong Zhang 
447b5b4017aSHong Zhang   Level: developer
448b5b4017aSHong Zhang 
449bcf0153eSBarry Smith   Note:
450bcf0153eSBarry Smith   `TSComputeIHessianProductFunctionUU()` is typically used for sensitivity implementation,
451bcf0153eSBarry Smith   so most users would not generally call this routine themselves.
452bcf0153eSBarry Smith 
4531cc06b55SBarry Smith .seealso: [](ch_ts), `TSSetIHessianProduct()`
454b5b4017aSHong Zhang @*/
455d71ae5a4SJacob Faibussowitsch PetscErrorCode TSComputeIHessianProductFunctionUU(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV)
456d71ae5a4SJacob Faibussowitsch {
457b5b4017aSHong Zhang   PetscFunctionBegin;
4583ba16761SJacob Faibussowitsch   if (!VHV) PetscFunctionReturn(PETSC_SUCCESS);
459b5b4017aSHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
460b5b4017aSHong Zhang   PetscValidHeaderSpecific(U, VEC_CLASSID, 3);
461b5b4017aSHong Zhang 
462792fecdfSBarry Smith   if (ts->ihessianproduct_fuu) PetscCallBack("TS callback IHessianProduct 1 for sensitivity analysis", (*ts->ihessianproduct_fuu)(ts, t, U, Vl, Vr, VHV, ts->ihessianproductctx));
463ef1023bdSBarry Smith 
46467633408SHong Zhang   /* does not consider IMEX for now, so either IHessian or RHSHessian will be calculated, using the same output VHV */
46567633408SHong Zhang   if (ts->rhshessianproduct_guu) {
46667633408SHong Zhang     PetscInt nadj;
4679566063dSJacob Faibussowitsch     PetscCall(TSComputeRHSHessianProductFunctionUU(ts, t, U, Vl, Vr, VHV));
46848a46eb9SPierre Jolivet     for (nadj = 0; nadj < ts->numcost; nadj++) PetscCall(VecScale(VHV[nadj], -1));
46967633408SHong Zhang   }
4703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
471b5b4017aSHong Zhang }
472b5b4017aSHong Zhang 
473b5b4017aSHong Zhang /*@C
474b1e111ebSHong Zhang   TSComputeIHessianProductFunctionUP - Runs the user-defined vector-Hessian-vector product function for Fup.
475b5b4017aSHong Zhang 
476c3339decSBarry Smith   Collective
477b5b4017aSHong Zhang 
478b5b4017aSHong Zhang   Input Parameters:
4792fe279fdSBarry Smith + ts - The `TS` context obtained from `TSCreate()`
4802fe279fdSBarry Smith . t  - the time
4812fe279fdSBarry Smith . U  - the solution at which to compute the Hessian product
4822fe279fdSBarry Smith . Vl - the array of input vectors to be multiplied with the Hessian from the left
4832fe279fdSBarry Smith - Vr - the input vector to be multiplied with the Hessian from the right
4842fe279fdSBarry Smith 
4852fe279fdSBarry Smith   Output Parameter:
486b43aa488SJacob Faibussowitsch . VHV - the array of output vectors that store the Hessian product
487b5b4017aSHong Zhang 
488b5b4017aSHong Zhang   Level: developer
489b5b4017aSHong Zhang 
490bcf0153eSBarry Smith   Note:
491bcf0153eSBarry Smith   `TSComputeIHessianProductFunctionUP()` is typically used for sensitivity implementation,
492bcf0153eSBarry Smith   so most users would not generally call this routine themselves.
493bcf0153eSBarry Smith 
4941cc06b55SBarry Smith .seealso: [](ch_ts), `TSSetIHessianProduct()`
495b5b4017aSHong Zhang @*/
496d71ae5a4SJacob Faibussowitsch PetscErrorCode TSComputeIHessianProductFunctionUP(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV)
497d71ae5a4SJacob Faibussowitsch {
498b5b4017aSHong Zhang   PetscFunctionBegin;
4993ba16761SJacob Faibussowitsch   if (!VHV) PetscFunctionReturn(PETSC_SUCCESS);
500b5b4017aSHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
501b5b4017aSHong Zhang   PetscValidHeaderSpecific(U, VEC_CLASSID, 3);
502b5b4017aSHong Zhang 
503792fecdfSBarry Smith   if (ts->ihessianproduct_fup) PetscCallBack("TS callback IHessianProduct 2 for sensitivity analysis", (*ts->ihessianproduct_fup)(ts, t, U, Vl, Vr, VHV, ts->ihessianproductctx));
504ef1023bdSBarry Smith 
50567633408SHong Zhang   /* does not consider IMEX for now, so either IHessian or RHSHessian will be calculated, using the same output VHV */
50667633408SHong Zhang   if (ts->rhshessianproduct_gup) {
50767633408SHong Zhang     PetscInt nadj;
5089566063dSJacob Faibussowitsch     PetscCall(TSComputeRHSHessianProductFunctionUP(ts, t, U, Vl, Vr, VHV));
50948a46eb9SPierre Jolivet     for (nadj = 0; nadj < ts->numcost; nadj++) PetscCall(VecScale(VHV[nadj], -1));
51067633408SHong Zhang   }
5113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
512b5b4017aSHong Zhang }
513b5b4017aSHong Zhang 
514b5b4017aSHong Zhang /*@C
515b1e111ebSHong Zhang   TSComputeIHessianProductFunctionPU - Runs the user-defined vector-Hessian-vector product function for Fpu.
516b5b4017aSHong Zhang 
517c3339decSBarry Smith   Collective
518b5b4017aSHong Zhang 
519b5b4017aSHong Zhang   Input Parameters:
5202fe279fdSBarry Smith + ts - The `TS` context obtained from `TSCreate()`
5212fe279fdSBarry Smith . t  - the time
5222fe279fdSBarry Smith . U  - the solution at which to compute the Hessian product
5232fe279fdSBarry Smith . Vl - the array of input vectors to be multiplied with the Hessian from the left
5242fe279fdSBarry Smith - Vr - the input vector to be multiplied with the Hessian from the right
5252fe279fdSBarry Smith 
5262fe279fdSBarry Smith   Output Parameter:
527b43aa488SJacob Faibussowitsch . VHV - the array of output vectors that store the Hessian product
528b5b4017aSHong Zhang 
529b5b4017aSHong Zhang   Level: developer
530b5b4017aSHong Zhang 
531bcf0153eSBarry Smith   Note:
532bcf0153eSBarry Smith   `TSComputeIHessianProductFunctionPU()` is typically used for sensitivity implementation,
533bcf0153eSBarry Smith   so most users would not generally call this routine themselves.
534bcf0153eSBarry Smith 
5351cc06b55SBarry Smith .seealso: [](ch_ts), `TSSetIHessianProduct()`
536b5b4017aSHong Zhang @*/
537d71ae5a4SJacob Faibussowitsch PetscErrorCode TSComputeIHessianProductFunctionPU(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV)
538d71ae5a4SJacob Faibussowitsch {
539b5b4017aSHong Zhang   PetscFunctionBegin;
5403ba16761SJacob Faibussowitsch   if (!VHV) PetscFunctionReturn(PETSC_SUCCESS);
541b5b4017aSHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
542b5b4017aSHong Zhang   PetscValidHeaderSpecific(U, VEC_CLASSID, 3);
543b5b4017aSHong Zhang 
544792fecdfSBarry Smith   if (ts->ihessianproduct_fpu) PetscCallBack("TS callback IHessianProduct 3 for sensitivity analysis", (*ts->ihessianproduct_fpu)(ts, t, U, Vl, Vr, VHV, ts->ihessianproductctx));
545ef1023bdSBarry Smith 
54667633408SHong Zhang   /* does not consider IMEX for now, so either IHessian or RHSHessian will be calculated, using the same output VHV */
54767633408SHong Zhang   if (ts->rhshessianproduct_gpu) {
54867633408SHong Zhang     PetscInt nadj;
5499566063dSJacob Faibussowitsch     PetscCall(TSComputeRHSHessianProductFunctionPU(ts, t, U, Vl, Vr, VHV));
55048a46eb9SPierre Jolivet     for (nadj = 0; nadj < ts->numcost; nadj++) PetscCall(VecScale(VHV[nadj], -1));
55167633408SHong Zhang   }
5523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
553b5b4017aSHong Zhang }
554b5b4017aSHong Zhang 
555b5b4017aSHong Zhang /*@C
556b1e111ebSHong Zhang   TSComputeIHessianProductFunctionPP - Runs the user-defined vector-Hessian-vector product function for Fpp.
557b5b4017aSHong Zhang 
558c3339decSBarry Smith   Collective
559b5b4017aSHong Zhang 
560b5b4017aSHong Zhang   Input Parameters:
5612fe279fdSBarry Smith + ts - The `TS` context obtained from `TSCreate()`
5622fe279fdSBarry Smith . t  - the time
5632fe279fdSBarry Smith . U  - the solution at which to compute the Hessian product
5642fe279fdSBarry Smith . Vl - the array of input vectors to be multiplied with the Hessian from the left
5652fe279fdSBarry Smith - Vr - the input vector to be multiplied with the Hessian from the right
5662fe279fdSBarry Smith 
5672fe279fdSBarry Smith   Output Parameter:
568b43aa488SJacob Faibussowitsch . VHV - the array of output vectors that store the Hessian product
569b5b4017aSHong Zhang 
570b5b4017aSHong Zhang   Level: developer
571b5b4017aSHong Zhang 
572bcf0153eSBarry Smith   Note:
573bcf0153eSBarry Smith   `TSComputeIHessianProductFunctionPP()` is typically used for sensitivity implementation,
574bcf0153eSBarry Smith   so most users would not generally call this routine themselves.
575bcf0153eSBarry Smith 
5761cc06b55SBarry Smith .seealso: [](ch_ts), `TSSetIHessianProduct()`
577b5b4017aSHong Zhang @*/
578d71ae5a4SJacob Faibussowitsch PetscErrorCode TSComputeIHessianProductFunctionPP(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV)
579d71ae5a4SJacob Faibussowitsch {
580b5b4017aSHong Zhang   PetscFunctionBegin;
5813ba16761SJacob Faibussowitsch   if (!VHV) PetscFunctionReturn(PETSC_SUCCESS);
582b5b4017aSHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
583b5b4017aSHong Zhang   PetscValidHeaderSpecific(U, VEC_CLASSID, 3);
584b5b4017aSHong Zhang 
585792fecdfSBarry Smith   if (ts->ihessianproduct_fpp) PetscCallBack("TS callback IHessianProduct 3 for sensitivity analysis", (*ts->ihessianproduct_fpp)(ts, t, U, Vl, Vr, VHV, ts->ihessianproductctx));
586ef1023bdSBarry Smith 
58767633408SHong Zhang   /* does not consider IMEX for now, so either IHessian or RHSHessian will be calculated, using the same output VHV */
58867633408SHong Zhang   if (ts->rhshessianproduct_gpp) {
58967633408SHong Zhang     PetscInt nadj;
5909566063dSJacob Faibussowitsch     PetscCall(TSComputeRHSHessianProductFunctionPP(ts, t, U, Vl, Vr, VHV));
59148a46eb9SPierre Jolivet     for (nadj = 0; nadj < ts->numcost; nadj++) PetscCall(VecScale(VHV[nadj], -1));
59267633408SHong Zhang   }
5933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
594b5b4017aSHong Zhang }
595b5b4017aSHong Zhang 
59614d0ab18SJacob Faibussowitsch // PetscClangLinter pragma disable: -fdoc-param-list-func-parameter-documentation
59714d0ab18SJacob Faibussowitsch // PetscClangLinter pragma disable: -fdoc-section-header-unknown
59813af1a74SHong Zhang /*@C
59914d0ab18SJacob Faibussowitsch   TSSetRHSHessianProduct - Sets the function that computes the vector-Hessian-vector
60014d0ab18SJacob Faibussowitsch   product. The Hessian is the second-order derivative of G (RHSFunction) w.r.t. the state
60114d0ab18SJacob Faibussowitsch   variable.
60213af1a74SHong Zhang 
603c3339decSBarry Smith   Logically Collective
60413af1a74SHong Zhang 
60513af1a74SHong Zhang   Input Parameters:
606bcf0153eSBarry Smith + ts     - `TS` context obtained from `TSCreate()`
60713af1a74SHong Zhang . rhshp1 - an array of vectors storing the result of vector-Hessian-vector product for G_UU
60813af1a74SHong Zhang . hessianproductfunc1 - vector-Hessian-vector product function for G_UU
60913af1a74SHong Zhang . rhshp2 - an array of vectors storing the result of vector-Hessian-vector product for G_UP
61013af1a74SHong Zhang . hessianproductfunc2 - vector-Hessian-vector product function for G_UP
61113af1a74SHong Zhang . rhshp3 - an array of vectors storing the result of vector-Hessian-vector product for G_PU
61213af1a74SHong Zhang . hessianproductfunc3 - vector-Hessian-vector product function for G_PU
61313af1a74SHong Zhang . rhshp4 - an array of vectors storing the result of vector-Hessian-vector product for G_PP
61414d0ab18SJacob Faibussowitsch . hessianproductfunc4 - vector-Hessian-vector product function for G_PP
61514d0ab18SJacob Faibussowitsch - ctx    - [optional] user-defined function context
61613af1a74SHong Zhang 
61714d0ab18SJacob Faibussowitsch   Calling sequence of `rhshessianproductfunc1`:
61814d0ab18SJacob Faibussowitsch + ts  - the `TS` context
61914d0ab18SJacob Faibussowitsch . t   - current timestep
62013af1a74SHong Zhang . U   - input vector (current ODE solution)
621369cf35fSHong Zhang . Vl  - an array of input vectors to be left-multiplied with the Hessian
62213af1a74SHong Zhang . Vr  - input vector to be right-multiplied with the Hessian
623369cf35fSHong Zhang . VHV - an array of output vectors for vector-Hessian-vector product
62413af1a74SHong Zhang - ctx - [optional] user-defined function context
62513af1a74SHong Zhang 
62613af1a74SHong Zhang   Level: intermediate
62713af1a74SHong Zhang 
628369cf35fSHong Zhang   Notes:
62914d0ab18SJacob Faibussowitsch   All other functions have the same calling sequence as `rhhessianproductfunc1`, so their
63014d0ab18SJacob Faibussowitsch   descriptions are omitted for brevity.
63114d0ab18SJacob Faibussowitsch 
632369cf35fSHong Zhang   The first Hessian function and the working array are required.
63314d0ab18SJacob Faibussowitsch 
634369cf35fSHong Zhang   As an example to implement the callback functions, the second callback function calculates the vector-Hessian-vector product
635369cf35fSHong Zhang   $ Vl_n^T*G_UP*Vr
636369cf35fSHong 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.
637369cf35fSHong Zhang   Each entry of G_UP corresponds to the derivative
638369cf35fSHong Zhang   $ G_UP[i][j][k] = \frac{\partial^2 G[i]}{\partial U[j] \partial P[k]}.
639369cf35fSHong 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
640369cf35fSHong Zhang   $ VHV_n[j] = \sum_i \sum_k {Vl_n[i] * G_UP[i][j][k] * Vr[k]}
641369cf35fSHong Zhang   If the cost function is a scalar, there will be only one vector in Vl and VHV.
64213af1a74SHong Zhang 
6432fe279fdSBarry Smith .seealso: `TS`, `TSAdjoint`
64413af1a74SHong Zhang @*/
64514d0ab18SJacob Faibussowitsch PetscErrorCode TSSetRHSHessianProduct(TS ts, Vec *rhshp1, PetscErrorCode (*rhshessianproductfunc1)(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV, void *ctx), 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)
646d71ae5a4SJacob Faibussowitsch {
64713af1a74SHong Zhang   PetscFunctionBegin;
64813af1a74SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
6494f572ea9SToby Isaac   PetscAssertPointer(rhshp1, 2);
65013af1a74SHong Zhang 
65113af1a74SHong Zhang   ts->rhshessianproductctx = ctx;
65213af1a74SHong Zhang   if (rhshp1) ts->vecs_guu = rhshp1;
65313af1a74SHong Zhang   if (rhshp2) ts->vecs_gup = rhshp2;
65413af1a74SHong Zhang   if (rhshp3) ts->vecs_gpu = rhshp3;
65513af1a74SHong Zhang   if (rhshp4) ts->vecs_gpp = rhshp4;
65613af1a74SHong Zhang   ts->rhshessianproduct_guu = rhshessianproductfunc1;
65713af1a74SHong Zhang   ts->rhshessianproduct_gup = rhshessianproductfunc2;
65813af1a74SHong Zhang   ts->rhshessianproduct_gpu = rhshessianproductfunc3;
65913af1a74SHong Zhang   ts->rhshessianproduct_gpp = rhshessianproductfunc4;
6603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
66113af1a74SHong Zhang }
66213af1a74SHong Zhang 
66313af1a74SHong Zhang /*@C
664b1e111ebSHong Zhang   TSComputeRHSHessianProductFunctionUU - Runs the user-defined vector-Hessian-vector product function for Guu.
66513af1a74SHong Zhang 
666c3339decSBarry Smith   Collective
66713af1a74SHong Zhang 
66813af1a74SHong Zhang   Input Parameters:
6692fe279fdSBarry Smith + ts - The `TS` context obtained from `TSCreate()`
6702fe279fdSBarry Smith . t  - the time
6712fe279fdSBarry Smith . U  - the solution at which to compute the Hessian product
6722fe279fdSBarry Smith . Vl - the array of input vectors to be multiplied with the Hessian from the left
6732fe279fdSBarry Smith - Vr - the input vector to be multiplied with the Hessian from the right
6742fe279fdSBarry Smith 
6752fe279fdSBarry Smith   Output Parameter:
676b43aa488SJacob Faibussowitsch . VHV - the array of output vectors that store the Hessian product
67713af1a74SHong Zhang 
67813af1a74SHong Zhang   Level: developer
67913af1a74SHong Zhang 
680bcf0153eSBarry Smith   Note:
681bcf0153eSBarry Smith   `TSComputeRHSHessianProductFunctionUU()` is typically used for sensitivity implementation,
682bcf0153eSBarry Smith   so most users would not generally call this routine themselves.
683bcf0153eSBarry Smith 
6841cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSetRHSHessianProduct()`
68513af1a74SHong Zhang @*/
686d71ae5a4SJacob Faibussowitsch PetscErrorCode TSComputeRHSHessianProductFunctionUU(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV)
687d71ae5a4SJacob Faibussowitsch {
68813af1a74SHong Zhang   PetscFunctionBegin;
6893ba16761SJacob Faibussowitsch   if (!VHV) PetscFunctionReturn(PETSC_SUCCESS);
69013af1a74SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
69113af1a74SHong Zhang   PetscValidHeaderSpecific(U, VEC_CLASSID, 3);
69213af1a74SHong Zhang 
693792fecdfSBarry Smith   PetscCallBack("TS callback RHSHessianProduct 1 for sensitivity analysis", (*ts->rhshessianproduct_guu)(ts, t, U, Vl, Vr, VHV, ts->rhshessianproductctx));
6943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
69513af1a74SHong Zhang }
69613af1a74SHong Zhang 
69713af1a74SHong Zhang /*@C
698b1e111ebSHong Zhang   TSComputeRHSHessianProductFunctionUP - Runs the user-defined vector-Hessian-vector product function for Gup.
69913af1a74SHong Zhang 
700c3339decSBarry Smith   Collective
70113af1a74SHong Zhang 
70213af1a74SHong Zhang   Input Parameters:
7032fe279fdSBarry Smith + ts - The `TS` context obtained from `TSCreate()`
7042fe279fdSBarry Smith . t  - the time
7052fe279fdSBarry Smith . U  - the solution at which to compute the Hessian product
7062fe279fdSBarry Smith . Vl - the array of input vectors to be multiplied with the Hessian from the left
7072fe279fdSBarry Smith - Vr - the input vector to be multiplied with the Hessian from the right
7082fe279fdSBarry Smith 
7092fe279fdSBarry Smith   Output Parameter:
710b43aa488SJacob Faibussowitsch . VHV - the array of output vectors that store the Hessian product
71113af1a74SHong Zhang 
71213af1a74SHong Zhang   Level: developer
71313af1a74SHong Zhang 
714bcf0153eSBarry Smith   Note:
715bcf0153eSBarry Smith   `TSComputeRHSHessianProductFunctionUP()` is typically used for sensitivity implementation,
716bcf0153eSBarry Smith   so most users would not generally call this routine themselves.
717bcf0153eSBarry Smith 
7181cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSetRHSHessianProduct()`
71913af1a74SHong Zhang @*/
720d71ae5a4SJacob Faibussowitsch PetscErrorCode TSComputeRHSHessianProductFunctionUP(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV)
721d71ae5a4SJacob Faibussowitsch {
72213af1a74SHong Zhang   PetscFunctionBegin;
7233ba16761SJacob Faibussowitsch   if (!VHV) PetscFunctionReturn(PETSC_SUCCESS);
72413af1a74SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
72513af1a74SHong Zhang   PetscValidHeaderSpecific(U, VEC_CLASSID, 3);
72613af1a74SHong Zhang 
727792fecdfSBarry Smith   PetscCallBack("TS callback RHSHessianProduct 2 for sensitivity analysis", (*ts->rhshessianproduct_gup)(ts, t, U, Vl, Vr, VHV, ts->rhshessianproductctx));
7283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
72913af1a74SHong Zhang }
73013af1a74SHong Zhang 
73113af1a74SHong Zhang /*@C
732b1e111ebSHong Zhang   TSComputeRHSHessianProductFunctionPU - Runs the user-defined vector-Hessian-vector product function for Gpu.
73313af1a74SHong Zhang 
734c3339decSBarry Smith   Collective
73513af1a74SHong Zhang 
73613af1a74SHong Zhang   Input Parameters:
7372fe279fdSBarry Smith + ts - The `TS` context obtained from `TSCreate()`
7382fe279fdSBarry Smith . t  - the time
7392fe279fdSBarry Smith . U  - the solution at which to compute the Hessian product
7402fe279fdSBarry Smith . Vl - the array of input vectors to be multiplied with the Hessian from the left
7412fe279fdSBarry Smith - Vr - the input vector to be multiplied with the Hessian from the right
7422fe279fdSBarry Smith 
7432fe279fdSBarry Smith   Output Parameter:
744b43aa488SJacob Faibussowitsch . VHV - the array of output vectors that store the Hessian product
74513af1a74SHong Zhang 
74613af1a74SHong Zhang   Level: developer
74713af1a74SHong Zhang 
748bcf0153eSBarry Smith   Note:
749bcf0153eSBarry Smith   `TSComputeRHSHessianProductFunctionPU()` is typically used for sensitivity implementation,
750bcf0153eSBarry Smith   so most users would not generally call this routine themselves.
751bcf0153eSBarry Smith 
7521cc06b55SBarry Smith .seealso: [](ch_ts), `TSSetRHSHessianProduct()`
75313af1a74SHong Zhang @*/
754d71ae5a4SJacob Faibussowitsch PetscErrorCode TSComputeRHSHessianProductFunctionPU(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV)
755d71ae5a4SJacob Faibussowitsch {
75613af1a74SHong Zhang   PetscFunctionBegin;
7573ba16761SJacob Faibussowitsch   if (!VHV) PetscFunctionReturn(PETSC_SUCCESS);
75813af1a74SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
75913af1a74SHong Zhang   PetscValidHeaderSpecific(U, VEC_CLASSID, 3);
76013af1a74SHong Zhang 
761792fecdfSBarry Smith   PetscCallBack("TS callback RHSHessianProduct 3 for sensitivity analysis", (*ts->rhshessianproduct_gpu)(ts, t, U, Vl, Vr, VHV, ts->rhshessianproductctx));
7623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
76313af1a74SHong Zhang }
76413af1a74SHong Zhang 
76513af1a74SHong Zhang /*@C
766b1e111ebSHong Zhang   TSComputeRHSHessianProductFunctionPP - Runs the user-defined vector-Hessian-vector product function for Gpp.
76713af1a74SHong Zhang 
768c3339decSBarry Smith   Collective
76913af1a74SHong Zhang 
77013af1a74SHong Zhang   Input Parameters:
7712fe279fdSBarry Smith + ts - The `TS` context obtained from `TSCreate()`
7722fe279fdSBarry Smith . t  - the time
7732fe279fdSBarry Smith . U  - the solution at which to compute the Hessian product
7742fe279fdSBarry Smith . Vl - the array of input vectors to be multiplied with the Hessian from the left
7752fe279fdSBarry Smith - Vr - the input vector to be multiplied with the Hessian from the right
7762fe279fdSBarry Smith 
7772fe279fdSBarry Smith   Output Parameter:
778b43aa488SJacob Faibussowitsch . VHV - the array of output vectors that store the Hessian product
77913af1a74SHong Zhang 
78013af1a74SHong Zhang   Level: developer
78113af1a74SHong Zhang 
782bcf0153eSBarry Smith   Note:
783bcf0153eSBarry Smith   `TSComputeRHSHessianProductFunctionPP()` is typically used for sensitivity implementation,
784bcf0153eSBarry Smith   so most users would not generally call this routine themselves.
785bcf0153eSBarry Smith 
7861cc06b55SBarry Smith .seealso: [](ch_ts), `TSSetRHSHessianProduct()`
78713af1a74SHong Zhang @*/
788d71ae5a4SJacob Faibussowitsch PetscErrorCode TSComputeRHSHessianProductFunctionPP(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV)
789d71ae5a4SJacob Faibussowitsch {
79013af1a74SHong Zhang   PetscFunctionBegin;
7913ba16761SJacob Faibussowitsch   if (!VHV) PetscFunctionReturn(PETSC_SUCCESS);
79213af1a74SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
79313af1a74SHong Zhang   PetscValidHeaderSpecific(U, VEC_CLASSID, 3);
79413af1a74SHong Zhang 
795792fecdfSBarry Smith   PetscCallBack("TS callback RHSHessianProduct 3 for sensitivity analysis", (*ts->rhshessianproduct_gpp)(ts, t, U, Vl, Vr, VHV, ts->rhshessianproductctx));
7963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
79713af1a74SHong Zhang }
79813af1a74SHong Zhang 
799814e85d6SHong Zhang /* --------------------------- Adjoint sensitivity ---------------------------*/
800814e85d6SHong Zhang 
801814e85d6SHong Zhang /*@
802814e85d6SHong 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
803bcf0153eSBarry Smith   for use by the `TS` adjoint routines.
804814e85d6SHong Zhang 
805c3339decSBarry Smith   Logically Collective
806814e85d6SHong Zhang 
807814e85d6SHong Zhang   Input Parameters:
808bcf0153eSBarry Smith + ts      - the `TS` context obtained from `TSCreate()`
809d2fd7bfcSHong Zhang . numcost - number of gradients to be computed, this is the number of cost functions
810814e85d6SHong 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
811814e85d6SHong Zhang - mu      - gradients with respect to the parameters, the number of entries in these vectors is the same as the number of parameters
812814e85d6SHong Zhang 
813814e85d6SHong Zhang   Level: beginner
814814e85d6SHong Zhang 
815814e85d6SHong Zhang   Notes:
81635cb6cd3SPierre Jolivet   the entries in these vectors must be correctly initialized with the values lambda_i = df/dy|finaltime  mu_i = df/dp|finaltime
817814e85d6SHong Zhang 
81835cb6cd3SPierre Jolivet   After `TSAdjointSolve()` is called the lambda and the mu contain the computed sensitivities
819814e85d6SHong Zhang 
820bcf0153eSBarry Smith .seealso: `TS`, `TSAdjointSolve()`, `TSGetCostGradients()`
821814e85d6SHong Zhang @*/
822d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSetCostGradients(TS ts, PetscInt numcost, Vec *lambda, Vec *mu)
823d71ae5a4SJacob Faibussowitsch {
824814e85d6SHong Zhang   PetscFunctionBegin;
825814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
8264f572ea9SToby Isaac   PetscAssertPointer(lambda, 3);
827814e85d6SHong Zhang   ts->vecs_sensi  = lambda;
828814e85d6SHong Zhang   ts->vecs_sensip = mu;
8293c633725SBarry 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");
830814e85d6SHong Zhang   ts->numcost = numcost;
8313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
832814e85d6SHong Zhang }
833814e85d6SHong Zhang 
834814e85d6SHong Zhang /*@
835bcf0153eSBarry Smith   TSGetCostGradients - Returns the gradients from the `TSAdjointSolve()`
836814e85d6SHong Zhang 
837bcf0153eSBarry Smith   Not Collective, but the vectors returned are parallel if `TS` is parallel
838814e85d6SHong Zhang 
839814e85d6SHong Zhang   Input Parameter:
840bcf0153eSBarry Smith . ts - the `TS` context obtained from `TSCreate()`
841814e85d6SHong Zhang 
842d8d19677SJose E. Roman   Output Parameters:
8436b867d5aSJose E. Roman + numcost - size of returned arrays
8446b867d5aSJose E. Roman . lambda  - vectors containing the gradients of the cost functions with respect to the ODE/DAE solution variables
845814e85d6SHong Zhang - mu      - vectors containing the gradients of the cost functions with respect to the problem parameters
846814e85d6SHong Zhang 
847814e85d6SHong Zhang   Level: intermediate
848814e85d6SHong Zhang 
8491cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSAdjointSolve()`, `TSSetCostGradients()`
850814e85d6SHong Zhang @*/
851d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGetCostGradients(TS ts, PetscInt *numcost, Vec **lambda, Vec **mu)
852d71ae5a4SJacob Faibussowitsch {
853814e85d6SHong Zhang   PetscFunctionBegin;
854814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
855814e85d6SHong Zhang   if (numcost) *numcost = ts->numcost;
856814e85d6SHong Zhang   if (lambda) *lambda = ts->vecs_sensi;
857814e85d6SHong Zhang   if (mu) *mu = ts->vecs_sensip;
8583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
859814e85d6SHong Zhang }
860814e85d6SHong Zhang 
861814e85d6SHong Zhang /*@
862b5b4017aSHong 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
863bcf0153eSBarry Smith   for use by the `TS` adjoint routines.
864b5b4017aSHong Zhang 
865c3339decSBarry Smith   Logically Collective
866b5b4017aSHong Zhang 
867b5b4017aSHong Zhang   Input Parameters:
868bcf0153eSBarry Smith + ts      - the `TS` context obtained from `TSCreate()`
869b5b4017aSHong Zhang . numcost - number of cost functions
870b5b4017aSHong 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
871b5b4017aSHong 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
872b5b4017aSHong Zhang - dir     - the direction vector that are multiplied with the Hessian of the cost functions
873b5b4017aSHong Zhang 
874b5b4017aSHong Zhang   Level: beginner
875b5b4017aSHong Zhang 
876bcf0153eSBarry Smith   Notes:
877bcf0153eSBarry Smith   Hessian of the cost function is completely different from Hessian of the ODE/DAE system
878b5b4017aSHong Zhang 
879bcf0153eSBarry Smith   For second-order adjoint, one needs to call this function and then `TSAdjointSetForward()` before `TSSolve()`.
880b5b4017aSHong Zhang 
88135cb6cd3SPierre 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.
882b5b4017aSHong Zhang 
8832fe279fdSBarry Smith   Passing `NULL` for `lambda2` disables the second-order calculation.
884bcf0153eSBarry Smith 
8851cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSAdjointSolve()`, `TSAdjointSetForward()`
886b5b4017aSHong Zhang @*/
887d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSetCostHessianProducts(TS ts, PetscInt numcost, Vec *lambda2, Vec *mu2, Vec dir)
888d71ae5a4SJacob Faibussowitsch {
889b5b4017aSHong Zhang   PetscFunctionBegin;
890b5b4017aSHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
8913c633725SBarry 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");
892b5b4017aSHong Zhang   ts->numcost      = numcost;
893b5b4017aSHong Zhang   ts->vecs_sensi2  = lambda2;
89433f72c5dSHong Zhang   ts->vecs_sensi2p = mu2;
895b5b4017aSHong Zhang   ts->vec_dir      = dir;
8963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
897b5b4017aSHong Zhang }
898b5b4017aSHong Zhang 
899b5b4017aSHong Zhang /*@
900bcf0153eSBarry Smith   TSGetCostHessianProducts - Returns the gradients from the `TSAdjointSolve()`
901b5b4017aSHong Zhang 
902bcf0153eSBarry Smith   Not Collective, but vectors returned are parallel if `TS` is parallel
903b5b4017aSHong Zhang 
904b5b4017aSHong Zhang   Input Parameter:
905bcf0153eSBarry Smith . ts - the `TS` context obtained from `TSCreate()`
906b5b4017aSHong Zhang 
907d8d19677SJose E. Roman   Output Parameters:
908b5b4017aSHong Zhang + numcost - number of cost functions
909b5b4017aSHong 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
910b5b4017aSHong 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
911b5b4017aSHong Zhang - dir     - the direction vector that are multiplied with the Hessian of the cost functions
912b5b4017aSHong Zhang 
913b5b4017aSHong Zhang   Level: intermediate
914b5b4017aSHong Zhang 
9151cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdjointSolve()`, `TSSetCostHessianProducts()`
916b5b4017aSHong Zhang @*/
917d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGetCostHessianProducts(TS ts, PetscInt *numcost, Vec **lambda2, Vec **mu2, Vec *dir)
918d71ae5a4SJacob Faibussowitsch {
919b5b4017aSHong Zhang   PetscFunctionBegin;
920b5b4017aSHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
921b5b4017aSHong Zhang   if (numcost) *numcost = ts->numcost;
922b5b4017aSHong Zhang   if (lambda2) *lambda2 = ts->vecs_sensi2;
92333f72c5dSHong Zhang   if (mu2) *mu2 = ts->vecs_sensi2p;
924b5b4017aSHong Zhang   if (dir) *dir = ts->vec_dir;
9253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
926b5b4017aSHong Zhang }
927b5b4017aSHong Zhang 
928b5b4017aSHong Zhang /*@
929ecf68647SHong Zhang   TSAdjointSetForward - Trigger the tangent linear solver and initialize the forward sensitivities
930b5b4017aSHong Zhang 
931c3339decSBarry Smith   Logically Collective
932b5b4017aSHong Zhang 
933b5b4017aSHong Zhang   Input Parameters:
934bcf0153eSBarry Smith + ts   - the `TS` context obtained from `TSCreate()`
935b5b4017aSHong Zhang - didp - the derivative of initial values w.r.t. parameters
936b5b4017aSHong Zhang 
937b5b4017aSHong Zhang   Level: intermediate
938b5b4017aSHong Zhang 
939bcf0153eSBarry Smith   Notes:
9402fe279fdSBarry 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.
9412fe279fdSBarry Smith   `TSAdjoint` does not reset the tangent linear solver automatically, `TSAdjointResetForward()` should be called to reset the tangent linear solver.
942b5b4017aSHong Zhang 
9431cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdjointSolve()`, `TSSetCostHessianProducts()`, `TSAdjointResetForward()`
944b5b4017aSHong Zhang @*/
945d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointSetForward(TS ts, Mat didp)
946d71ae5a4SJacob Faibussowitsch {
94733f72c5dSHong Zhang   Mat          A;
94833f72c5dSHong Zhang   Vec          sp;
94933f72c5dSHong Zhang   PetscScalar *xarr;
95033f72c5dSHong Zhang   PetscInt     lsize;
951b5b4017aSHong Zhang 
952b5b4017aSHong Zhang   PetscFunctionBegin;
953b5b4017aSHong Zhang   ts->forward_solve = PETSC_TRUE; /* turn on tangent linear mode */
9543c633725SBarry Smith   PetscCheck(ts->vecs_sensi2, PetscObjectComm((PetscObject)ts), PETSC_ERR_USER, "Must call TSSetCostHessianProducts() first");
9553c633725SBarry Smith   PetscCheck(ts->vec_dir, PetscObjectComm((PetscObject)ts), PETSC_ERR_USER, "Directional vector is missing. Call TSSetCostHessianProducts() to set it.");
95633f72c5dSHong Zhang   /* create a single-column dense matrix */
9579566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(ts->vec_sol, &lsize));
9589566063dSJacob Faibussowitsch   PetscCall(MatCreateDense(PetscObjectComm((PetscObject)ts), lsize, PETSC_DECIDE, PETSC_DECIDE, 1, NULL, &A));
95933f72c5dSHong Zhang 
9609566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &sp));
9619566063dSJacob Faibussowitsch   PetscCall(MatDenseGetColumn(A, 0, &xarr));
9629566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(sp, xarr));
963ecf68647SHong Zhang   if (ts->vecs_sensi2p) { /* tangent linear variable initialized as 2*dIdP*dir */
96433f72c5dSHong Zhang     if (didp) {
9659566063dSJacob Faibussowitsch       PetscCall(MatMult(didp, ts->vec_dir, sp));
9669566063dSJacob Faibussowitsch       PetscCall(VecScale(sp, 2.));
96733f72c5dSHong Zhang     } else {
9689566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(sp));
96933f72c5dSHong Zhang     }
970ecf68647SHong Zhang   } else { /* tangent linear variable initialized as dir */
9719566063dSJacob Faibussowitsch     PetscCall(VecCopy(ts->vec_dir, sp));
97233f72c5dSHong Zhang   }
9739566063dSJacob Faibussowitsch   PetscCall(VecResetArray(sp));
9749566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreColumn(A, &xarr));
9759566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&sp));
97633f72c5dSHong Zhang 
9779566063dSJacob Faibussowitsch   PetscCall(TSForwardSetInitialSensitivities(ts, A)); /* if didp is NULL, identity matrix is assumed */
97833f72c5dSHong Zhang 
9799566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
9803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
981b5b4017aSHong Zhang }
982b5b4017aSHong Zhang 
983b5b4017aSHong Zhang /*@
984ecf68647SHong Zhang   TSAdjointResetForward - Reset the tangent linear solver and destroy the tangent linear context
985ecf68647SHong Zhang 
986c3339decSBarry Smith   Logically Collective
987ecf68647SHong Zhang 
9882fe279fdSBarry Smith   Input Parameter:
989bcf0153eSBarry Smith . ts - the `TS` context obtained from `TSCreate()`
990ecf68647SHong Zhang 
991ecf68647SHong Zhang   Level: intermediate
992ecf68647SHong Zhang 
9931cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdjointSetForward()`
994ecf68647SHong Zhang @*/
995d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointResetForward(TS ts)
996d71ae5a4SJacob Faibussowitsch {
997ecf68647SHong Zhang   PetscFunctionBegin;
998ecf68647SHong Zhang   ts->forward_solve = PETSC_FALSE; /* turn off tangent linear mode */
9999566063dSJacob Faibussowitsch   PetscCall(TSForwardReset(ts));
10003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1001ecf68647SHong Zhang }
1002ecf68647SHong Zhang 
1003ecf68647SHong Zhang /*@
1004814e85d6SHong Zhang   TSAdjointSetUp - Sets up the internal data structures for the later use
1005814e85d6SHong Zhang   of an adjoint solver
1006814e85d6SHong Zhang 
1007c3339decSBarry Smith   Collective
1008814e85d6SHong Zhang 
1009814e85d6SHong Zhang   Input Parameter:
1010bcf0153eSBarry Smith . ts - the `TS` context obtained from `TSCreate()`
1011814e85d6SHong Zhang 
1012814e85d6SHong Zhang   Level: advanced
1013814e85d6SHong Zhang 
10141cc06b55SBarry Smith .seealso: [](ch_ts), `TSCreate()`, `TSAdjointStep()`, `TSSetCostGradients()`
1015814e85d6SHong Zhang @*/
1016d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointSetUp(TS ts)
1017d71ae5a4SJacob Faibussowitsch {
1018881c1a9bSHong Zhang   TSTrajectory tj;
1019881c1a9bSHong Zhang   PetscBool    match;
1020814e85d6SHong Zhang 
1021814e85d6SHong Zhang   PetscFunctionBegin;
1022814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
10233ba16761SJacob Faibussowitsch   if (ts->adjointsetupcalled) PetscFunctionReturn(PETSC_SUCCESS);
10243c633725SBarry Smith   PetscCheck(ts->vecs_sensi, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_WRONGSTATE, "Must call TSSetCostGradients() first");
10253c633725SBarry Smith   PetscCheck(!ts->vecs_sensip || ts->Jacp || ts->Jacprhs, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_WRONGSTATE, "Must call TSSetRHSJacobianP() or TSSetIJacobianP() first");
10269566063dSJacob Faibussowitsch   PetscCall(TSGetTrajectory(ts, &tj));
10279566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)tj, TSTRAJECTORYBASIC, &match));
1028881c1a9bSHong Zhang   if (match) {
1029881c1a9bSHong Zhang     PetscBool solution_only;
10309566063dSJacob Faibussowitsch     PetscCall(TSTrajectoryGetSolutionOnly(tj, &solution_only));
10313c633725SBarry 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");
1032881c1a9bSHong Zhang   }
10339566063dSJacob Faibussowitsch   PetscCall(TSTrajectorySetUseHistory(tj, PETSC_FALSE)); /* not use TSHistory */
103433f72c5dSHong Zhang 
1035cd4cee2dSHong Zhang   if (ts->quadraturets) { /* if there is integral in the cost function */
10369566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(ts->vecs_sensi[0], &ts->vec_drdu_col));
103748a46eb9SPierre Jolivet     if (ts->vecs_sensip) PetscCall(VecDuplicate(ts->vecs_sensip[0], &ts->vec_drdp_col));
1038814e85d6SHong Zhang   }
1039814e85d6SHong Zhang 
1040dbbe0bcdSBarry Smith   PetscTryTypeMethod(ts, adjointsetup);
1041814e85d6SHong Zhang   ts->adjointsetupcalled = PETSC_TRUE;
10423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1043814e85d6SHong Zhang }
1044814e85d6SHong Zhang 
1045814e85d6SHong Zhang /*@
1046bcf0153eSBarry Smith   TSAdjointReset - Resets a `TS` adjoint context and removes any allocated `Vec`s and `Mat`s.
1047ece46509SHong Zhang 
1048c3339decSBarry Smith   Collective
1049ece46509SHong Zhang 
1050ece46509SHong Zhang   Input Parameter:
1051bcf0153eSBarry Smith . ts - the `TS` context obtained from `TSCreate()`
1052ece46509SHong Zhang 
1053ece46509SHong Zhang   Level: beginner
1054ece46509SHong Zhang 
10551cc06b55SBarry Smith .seealso: [](ch_ts), `TSCreate()`, `TSAdjointSetUp()`, `TSADestroy()`
1056ece46509SHong Zhang @*/
1057d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointReset(TS ts)
1058d71ae5a4SJacob Faibussowitsch {
1059ece46509SHong Zhang   PetscFunctionBegin;
1060ece46509SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1061dbbe0bcdSBarry Smith   PetscTryTypeMethod(ts, adjointreset);
10627207e0fdSHong Zhang   if (ts->quadraturets) { /* if there is integral in the cost function */
10639566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&ts->vec_drdu_col));
106448a46eb9SPierre Jolivet     if (ts->vecs_sensip) PetscCall(VecDestroy(&ts->vec_drdp_col));
10657207e0fdSHong Zhang   }
1066ece46509SHong Zhang   ts->vecs_sensi         = NULL;
1067ece46509SHong Zhang   ts->vecs_sensip        = NULL;
1068ece46509SHong Zhang   ts->vecs_sensi2        = NULL;
106933f72c5dSHong Zhang   ts->vecs_sensi2p       = NULL;
1070ece46509SHong Zhang   ts->vec_dir            = NULL;
1071ece46509SHong Zhang   ts->adjointsetupcalled = PETSC_FALSE;
10723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1073ece46509SHong Zhang }
1074ece46509SHong Zhang 
1075ece46509SHong Zhang /*@
1076814e85d6SHong Zhang   TSAdjointSetSteps - Sets the number of steps the adjoint solver should take backward in time
1077814e85d6SHong Zhang 
1078c3339decSBarry Smith   Logically Collective
1079814e85d6SHong Zhang 
1080814e85d6SHong Zhang   Input Parameters:
1081bcf0153eSBarry Smith + ts    - the `TS` context obtained from `TSCreate()`
1082a2b725a8SWilliam Gropp - steps - number of steps to use
1083814e85d6SHong Zhang 
1084814e85d6SHong Zhang   Level: intermediate
1085814e85d6SHong Zhang 
1086814e85d6SHong Zhang   Notes:
1087bcf0153eSBarry Smith   Normally one does not call this and `TSAdjointSolve()` integrates back to the original timestep. One can call this
1088814e85d6SHong Zhang   so as to integrate back to less than the original timestep
1089814e85d6SHong Zhang 
10901cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdjointSolve()`, `TS`, `TSSetExactFinalTime()`
1091814e85d6SHong Zhang @*/
1092d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointSetSteps(TS ts, PetscInt steps)
1093d71ae5a4SJacob Faibussowitsch {
1094814e85d6SHong Zhang   PetscFunctionBegin;
1095814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1096814e85d6SHong Zhang   PetscValidLogicalCollectiveInt(ts, steps, 2);
10973c633725SBarry Smith   PetscCheck(steps >= 0, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_OUTOFRANGE, "Cannot step back a negative number of steps");
10983c633725SBarry Smith   PetscCheck(steps <= ts->steps, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_OUTOFRANGE, "Cannot step back more than the total number of forward steps");
1099814e85d6SHong Zhang   ts->adjoint_max_steps = steps;
11003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1101814e85d6SHong Zhang }
1102814e85d6SHong Zhang 
110314d0ab18SJacob Faibussowitsch // PetscClangLinter pragma disable: -fdoc-*
1104814e85d6SHong Zhang /*@C
1105bcf0153eSBarry Smith   TSAdjointSetRHSJacobian - Deprecated, use `TSSetRHSJacobianP()`
1106814e85d6SHong Zhang 
1107814e85d6SHong Zhang   Level: deprecated
1108814e85d6SHong Zhang @*/
1109d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointSetRHSJacobian(TS ts, Mat Amat, PetscErrorCode (*func)(TS, PetscReal, Vec, Mat, void *), void *ctx)
1110d71ae5a4SJacob Faibussowitsch {
1111814e85d6SHong Zhang   PetscFunctionBegin;
1112814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1113814e85d6SHong Zhang   PetscValidHeaderSpecific(Amat, MAT_CLASSID, 2);
1114814e85d6SHong Zhang 
1115814e85d6SHong Zhang   ts->rhsjacobianp    = func;
1116814e85d6SHong Zhang   ts->rhsjacobianpctx = ctx;
1117814e85d6SHong Zhang   if (Amat) {
11189566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)Amat));
11199566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&ts->Jacp));
1120814e85d6SHong Zhang     ts->Jacp = Amat;
1121814e85d6SHong Zhang   }
11223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1123814e85d6SHong Zhang }
1124814e85d6SHong Zhang 
112514d0ab18SJacob Faibussowitsch // PetscClangLinter pragma disable: -fdoc-*
1126814e85d6SHong Zhang /*@C
1127bcf0153eSBarry Smith   TSAdjointComputeRHSJacobian - Deprecated, use `TSComputeRHSJacobianP()`
1128814e85d6SHong Zhang 
1129814e85d6SHong Zhang   Level: deprecated
1130814e85d6SHong Zhang @*/
1131d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointComputeRHSJacobian(TS ts, PetscReal t, Vec U, Mat Amat)
1132d71ae5a4SJacob Faibussowitsch {
1133814e85d6SHong Zhang   PetscFunctionBegin;
1134814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1135c9ad9de0SHong Zhang   PetscValidHeaderSpecific(U, VEC_CLASSID, 3);
1136292bffcbSToby Isaac   PetscValidHeaderSpecific(Amat, MAT_CLASSID, 4);
1137814e85d6SHong Zhang 
1138792fecdfSBarry Smith   PetscCallBack("TS callback JacobianP for sensitivity analysis", (*ts->rhsjacobianp)(ts, t, U, Amat, ts->rhsjacobianpctx));
11393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1140814e85d6SHong Zhang }
1141814e85d6SHong Zhang 
114214d0ab18SJacob Faibussowitsch // PetscClangLinter pragma disable: -fdoc-*
1143814e85d6SHong Zhang /*@
1144bcf0153eSBarry Smith   TSAdjointComputeDRDYFunction - Deprecated, use `TSGetQuadratureTS()` then `TSComputeRHSJacobian()`
1145814e85d6SHong Zhang 
1146814e85d6SHong Zhang   Level: deprecated
1147814e85d6SHong Zhang @*/
1148d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointComputeDRDYFunction(TS ts, PetscReal t, Vec U, Vec *DRDU)
1149d71ae5a4SJacob Faibussowitsch {
1150814e85d6SHong Zhang   PetscFunctionBegin;
1151814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1152c9ad9de0SHong Zhang   PetscValidHeaderSpecific(U, VEC_CLASSID, 3);
1153814e85d6SHong Zhang 
1154792fecdfSBarry Smith   PetscCallBack("TS callback DRDY for sensitivity analysis", (*ts->drdufunction)(ts, t, U, DRDU, ts->costintegrandctx));
11553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1156814e85d6SHong Zhang }
1157814e85d6SHong Zhang 
115814d0ab18SJacob Faibussowitsch // PetscClangLinter pragma disable: -fdoc-*
1159814e85d6SHong Zhang /*@
1160bcf0153eSBarry Smith   TSAdjointComputeDRDPFunction - Deprecated, use `TSGetQuadratureTS()` then `TSComputeRHSJacobianP()`
1161814e85d6SHong Zhang 
1162814e85d6SHong Zhang   Level: deprecated
1163814e85d6SHong Zhang @*/
1164d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointComputeDRDPFunction(TS ts, PetscReal t, Vec U, Vec *DRDP)
1165d71ae5a4SJacob Faibussowitsch {
1166814e85d6SHong Zhang   PetscFunctionBegin;
1167814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1168c9ad9de0SHong Zhang   PetscValidHeaderSpecific(U, VEC_CLASSID, 3);
1169814e85d6SHong Zhang 
1170792fecdfSBarry Smith   PetscCallBack("TS callback DRDP for sensitivity analysis", (*ts->drdpfunction)(ts, t, U, DRDP, ts->costintegrandctx));
11713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1172814e85d6SHong Zhang }
1173814e85d6SHong Zhang 
117414d0ab18SJacob Faibussowitsch // PetscClangLinter pragma disable: -fdoc-param-list-func-parameter-documentation
1175814e85d6SHong Zhang /*@C
1176814e85d6SHong Zhang   TSAdjointMonitorSensi - monitors the first lambda sensitivity
1177814e85d6SHong Zhang 
1178814e85d6SHong Zhang   Level: intermediate
1179814e85d6SHong Zhang 
11801cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdjointMonitorSet()`
1181814e85d6SHong Zhang @*/
1182*ba38deedSJacob Faibussowitsch static PetscErrorCode TSAdjointMonitorSensi(TS ts, PetscInt step, PetscReal ptime, Vec v, PetscInt numcost, Vec *lambda, Vec *mu, PetscViewerAndFormat *vf)
1183d71ae5a4SJacob Faibussowitsch {
1184814e85d6SHong Zhang   PetscViewer viewer = vf->viewer;
1185814e85d6SHong Zhang 
1186814e85d6SHong Zhang   PetscFunctionBegin;
1187064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 8);
11889566063dSJacob Faibussowitsch   PetscCall(PetscViewerPushFormat(viewer, vf->format));
11899566063dSJacob Faibussowitsch   PetscCall(VecView(lambda[0], viewer));
11909566063dSJacob Faibussowitsch   PetscCall(PetscViewerPopFormat(viewer));
11913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1192814e85d6SHong Zhang }
1193814e85d6SHong Zhang 
1194814e85d6SHong Zhang /*@C
1195814e85d6SHong Zhang   TSAdjointMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user
1196814e85d6SHong Zhang 
1197c3339decSBarry Smith   Collective
1198814e85d6SHong Zhang 
1199814e85d6SHong Zhang   Input Parameters:
1200bcf0153eSBarry Smith + ts           - `TS` object you wish to monitor
1201814e85d6SHong Zhang . name         - the monitor type one is seeking
1202814e85d6SHong Zhang . help         - message indicating what monitoring is done
1203814e85d6SHong Zhang . manual       - manual page for the monitor
1204814e85d6SHong Zhang . monitor      - the monitor function
1205bcf0153eSBarry 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
1206814e85d6SHong Zhang 
1207814e85d6SHong Zhang   Level: developer
1208814e85d6SHong Zhang 
12091cc06b55SBarry Smith .seealso: [](ch_ts), `PetscOptionsGetViewer()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`,
1210db781477SPatrick Sanan           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`
1211b43aa488SJacob Faibussowitsch           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`,
1212db781477SPatrick Sanan           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
1213c2e3fba1SPatrick Sanan           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
1214db781477SPatrick Sanan           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
1215db781477SPatrick Sanan           `PetscOptionsFList()`, `PetscOptionsEList()`
1216814e85d6SHong Zhang @*/
1217d71ae5a4SJacob 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 *))
1218d71ae5a4SJacob Faibussowitsch {
1219814e85d6SHong Zhang   PetscViewer       viewer;
1220814e85d6SHong Zhang   PetscViewerFormat format;
1221814e85d6SHong Zhang   PetscBool         flg;
1222814e85d6SHong Zhang 
1223814e85d6SHong Zhang   PetscFunctionBegin;
12249566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)ts), ((PetscObject)ts)->options, ((PetscObject)ts)->prefix, name, &viewer, &format, &flg));
1225814e85d6SHong Zhang   if (flg) {
1226814e85d6SHong Zhang     PetscViewerAndFormat *vf;
12279566063dSJacob Faibussowitsch     PetscCall(PetscViewerAndFormatCreate(viewer, format, &vf));
12289566063dSJacob Faibussowitsch     PetscCall(PetscObjectDereference((PetscObject)viewer));
12291baa6e33SBarry Smith     if (monitorsetup) PetscCall((*monitorsetup)(ts, vf));
12309566063dSJacob Faibussowitsch     PetscCall(TSAdjointMonitorSet(ts, (PetscErrorCode(*)(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *, void *))monitor, vf, (PetscErrorCode(*)(void **))PetscViewerAndFormatDestroy));
1231814e85d6SHong Zhang   }
12323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1233814e85d6SHong Zhang }
1234814e85d6SHong Zhang 
1235814e85d6SHong Zhang /*@C
1236814e85d6SHong Zhang   TSAdjointMonitorSet - Sets an ADDITIONAL function that is to be used at every
1237814e85d6SHong Zhang   timestep to display the iteration's  progress.
1238814e85d6SHong Zhang 
1239c3339decSBarry Smith   Logically Collective
1240814e85d6SHong Zhang 
1241814e85d6SHong Zhang   Input Parameters:
1242bcf0153eSBarry Smith + ts              - the `TS` context obtained from `TSCreate()`
1243814e85d6SHong Zhang . adjointmonitor  - monitoring routine
124414d0ab18SJacob Faibussowitsch . adjointmctx     - [optional] user-defined context for private data for the monitor routine
124514d0ab18SJacob Faibussowitsch                     (use `NULL` if no context is desired)
124614d0ab18SJacob Faibussowitsch - adjointmdestroy - [optional] routine that frees monitor context (may be `NULL`)
1247814e85d6SHong Zhang 
124820f4b53cSBarry Smith   Calling sequence of `adjointmonitor`:
1249bcf0153eSBarry Smith + ts          - the `TS` context
125014d0ab18SJacob Faibussowitsch . steps       - iteration number (after the final time step the monitor routine is called with
125114d0ab18SJacob Faibussowitsch                 a step of -1, this is at the final time which may have been interpolated to)
1252814e85d6SHong Zhang . time        - current time
1253814e85d6SHong Zhang . u           - current iterate
1254814e85d6SHong Zhang . numcost     - number of cost functionos
1255814e85d6SHong Zhang . lambda      - sensitivities to initial conditions
1256814e85d6SHong Zhang . mu          - sensitivities to parameters
1257814e85d6SHong Zhang - adjointmctx - [optional] adjoint monitoring context
1258814e85d6SHong Zhang 
125914d0ab18SJacob Faibussowitsch   Calling sequence of `adjointmdestroy`:
126014d0ab18SJacob Faibussowitsch . mctx - the monitor context to destroy
126114d0ab18SJacob Faibussowitsch 
1262bcf0153eSBarry Smith   Level: intermediate
1263bcf0153eSBarry Smith 
1264bcf0153eSBarry Smith   Note:
1265814e85d6SHong Zhang   This routine adds an additional monitor to the list of monitors that
1266814e85d6SHong Zhang   already has been loaded.
1267814e85d6SHong Zhang 
1268b43aa488SJacob Faibussowitsch   Fortran Notes:
126920f4b53cSBarry Smith   Only a single monitor function can be set for each `TS` object
1270814e85d6SHong Zhang 
12711cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSAdjointSolve()`, `TSAdjointMonitorCancel()`
1272814e85d6SHong Zhang @*/
127314d0ab18SJacob Faibussowitsch PetscErrorCode TSAdjointMonitorSet(TS ts, PetscErrorCode (*adjointmonitor)(TS ts, PetscInt steps, PetscReal time, Vec u, PetscInt numcost, Vec *lambda, Vec *mu, void *adjointmctx), void *adjointmctx, PetscErrorCode (*adjointmdestroy)(void **mctx))
1274d71ae5a4SJacob Faibussowitsch {
1275814e85d6SHong Zhang   PetscInt  i;
1276814e85d6SHong Zhang   PetscBool identical;
1277814e85d6SHong Zhang 
1278814e85d6SHong Zhang   PetscFunctionBegin;
1279814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1280814e85d6SHong Zhang   for (i = 0; i < ts->numbermonitors; i++) {
12819566063dSJacob Faibussowitsch     PetscCall(PetscMonitorCompare((PetscErrorCode(*)(void))adjointmonitor, adjointmctx, adjointmdestroy, (PetscErrorCode(*)(void))ts->adjointmonitor[i], ts->adjointmonitorcontext[i], ts->adjointmonitordestroy[i], &identical));
12823ba16761SJacob Faibussowitsch     if (identical) PetscFunctionReturn(PETSC_SUCCESS);
1283814e85d6SHong Zhang   }
12843c633725SBarry Smith   PetscCheck(ts->numberadjointmonitors < MAXTSMONITORS, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many adjoint monitors set");
1285814e85d6SHong Zhang   ts->adjointmonitor[ts->numberadjointmonitors]          = adjointmonitor;
1286814e85d6SHong Zhang   ts->adjointmonitordestroy[ts->numberadjointmonitors]   = adjointmdestroy;
1287814e85d6SHong Zhang   ts->adjointmonitorcontext[ts->numberadjointmonitors++] = (void *)adjointmctx;
12883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1289814e85d6SHong Zhang }
1290814e85d6SHong Zhang 
1291814e85d6SHong Zhang /*@C
1292814e85d6SHong Zhang   TSAdjointMonitorCancel - Clears all the adjoint monitors that have been set on a time-step object.
1293814e85d6SHong Zhang 
1294c3339decSBarry Smith   Logically Collective
1295814e85d6SHong Zhang 
12962fe279fdSBarry Smith   Input Parameter:
1297bcf0153eSBarry Smith . ts - the `TS` context obtained from `TSCreate()`
1298814e85d6SHong Zhang 
1299814e85d6SHong Zhang   Notes:
1300814e85d6SHong Zhang   There is no way to remove a single, specific monitor.
1301814e85d6SHong Zhang 
1302814e85d6SHong Zhang   Level: intermediate
1303814e85d6SHong Zhang 
13041cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSAdjointSolve()`, `TSAdjointMonitorSet()`
1305814e85d6SHong Zhang @*/
1306d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointMonitorCancel(TS ts)
1307d71ae5a4SJacob Faibussowitsch {
1308814e85d6SHong Zhang   PetscInt i;
1309814e85d6SHong Zhang 
1310814e85d6SHong Zhang   PetscFunctionBegin;
1311814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1312814e85d6SHong Zhang   for (i = 0; i < ts->numberadjointmonitors; i++) {
131348a46eb9SPierre Jolivet     if (ts->adjointmonitordestroy[i]) PetscCall((*ts->adjointmonitordestroy[i])(&ts->adjointmonitorcontext[i]));
1314814e85d6SHong Zhang   }
1315814e85d6SHong Zhang   ts->numberadjointmonitors = 0;
13163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1317814e85d6SHong Zhang }
1318814e85d6SHong Zhang 
1319814e85d6SHong Zhang /*@C
1320814e85d6SHong Zhang   TSAdjointMonitorDefault - the default monitor of adjoint computations
1321814e85d6SHong Zhang 
132214d0ab18SJacob Faibussowitsch   Input Parameters:
132314d0ab18SJacob Faibussowitsch + ts      - the `TS` context
132414d0ab18SJacob Faibussowitsch . step    - iteration number (after the final time step the monitor routine is called with a
132514d0ab18SJacob Faibussowitsch step of -1, this is at the final time which may have been interpolated to)
132614d0ab18SJacob Faibussowitsch . time    - current time
132714d0ab18SJacob Faibussowitsch . v       - current iterate
132814d0ab18SJacob Faibussowitsch . numcost - number of cost functionos
132914d0ab18SJacob Faibussowitsch . lambda  - sensitivities to initial conditions
133014d0ab18SJacob Faibussowitsch . mu      - sensitivities to parameters
133114d0ab18SJacob Faibussowitsch - vf      - the viewer and format
133214d0ab18SJacob Faibussowitsch 
1333814e85d6SHong Zhang   Level: intermediate
1334814e85d6SHong Zhang 
13351cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSAdjointSolve()`, `TSAdjointMonitorSet()`
1336814e85d6SHong Zhang @*/
133714d0ab18SJacob Faibussowitsch PetscErrorCode TSAdjointMonitorDefault(TS ts, PetscInt step, PetscReal time, Vec v, PetscInt numcost, Vec *lambda, Vec *mu, PetscViewerAndFormat *vf)
1338d71ae5a4SJacob Faibussowitsch {
1339814e85d6SHong Zhang   PetscViewer viewer = vf->viewer;
1340814e85d6SHong Zhang 
1341814e85d6SHong Zhang   PetscFunctionBegin;
134214d0ab18SJacob Faibussowitsch   (void)v;
134314d0ab18SJacob Faibussowitsch   (void)numcost;
134414d0ab18SJacob Faibussowitsch   (void)lambda;
134514d0ab18SJacob Faibussowitsch   (void)mu;
1346064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 8);
13479566063dSJacob Faibussowitsch   PetscCall(PetscViewerPushFormat(viewer, vf->format));
13489566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)ts)->tablevel));
134914d0ab18SJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " TS dt %g time %g%s", step, (double)ts->time_step, (double)time, ts->steprollback ? " (r)\n" : "\n"));
13509566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)ts)->tablevel));
13519566063dSJacob Faibussowitsch   PetscCall(PetscViewerPopFormat(viewer));
13523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1353814e85d6SHong Zhang }
1354814e85d6SHong Zhang 
1355814e85d6SHong Zhang /*@C
1356bcf0153eSBarry Smith   TSAdjointMonitorDrawSensi - Monitors progress of the adjoint `TS` solvers by calling
1357bcf0153eSBarry Smith   `VecView()` for the sensitivities to initial states at each timestep
1358814e85d6SHong Zhang 
1359c3339decSBarry Smith   Collective
1360814e85d6SHong Zhang 
1361814e85d6SHong Zhang   Input Parameters:
1362bcf0153eSBarry Smith + ts      - the `TS` context
1363814e85d6SHong Zhang . step    - current time-step
1364814e85d6SHong Zhang . ptime   - current time
1365814e85d6SHong Zhang . u       - current state
1366814e85d6SHong Zhang . numcost - number of cost functions
1367814e85d6SHong Zhang . lambda  - sensitivities to initial conditions
1368814e85d6SHong Zhang . mu      - sensitivities to parameters
13692fe279fdSBarry Smith - dummy   - either a viewer or `NULL`
1370814e85d6SHong Zhang 
1371814e85d6SHong Zhang   Level: intermediate
1372814e85d6SHong Zhang 
13731cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdjointSolve()`, `TSAdjointMonitorSet()`, `TSAdjointMonitorDefault()`, `VecView()`
1374814e85d6SHong Zhang @*/
1375d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointMonitorDrawSensi(TS ts, PetscInt step, PetscReal ptime, Vec u, PetscInt numcost, Vec *lambda, Vec *mu, void *dummy)
1376d71ae5a4SJacob Faibussowitsch {
1377814e85d6SHong Zhang   TSMonitorDrawCtx ictx = (TSMonitorDrawCtx)dummy;
1378814e85d6SHong Zhang   PetscDraw        draw;
1379814e85d6SHong Zhang   PetscReal        xl, yl, xr, yr, h;
1380814e85d6SHong Zhang   char             time[32];
1381814e85d6SHong Zhang 
1382814e85d6SHong Zhang   PetscFunctionBegin;
13833ba16761SJacob Faibussowitsch   if (!(((ictx->howoften > 0) && (!(step % ictx->howoften))) || ((ictx->howoften == -1) && ts->reason))) PetscFunctionReturn(PETSC_SUCCESS);
1384814e85d6SHong Zhang 
13859566063dSJacob Faibussowitsch   PetscCall(VecView(lambda[0], ictx->viewer));
13869566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawGetDraw(ictx->viewer, 0, &draw));
13879566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(time, 32, "Timestep %d Time %g", (int)step, (double)ptime));
13889566063dSJacob Faibussowitsch   PetscCall(PetscDrawGetCoordinates(draw, &xl, &yl, &xr, &yr));
1389814e85d6SHong Zhang   h = yl + .95 * (yr - yl);
13909566063dSJacob Faibussowitsch   PetscCall(PetscDrawStringCentered(draw, .5 * (xl + xr), h, PETSC_DRAW_BLACK, time));
13919566063dSJacob Faibussowitsch   PetscCall(PetscDrawFlush(draw));
13923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1393814e85d6SHong Zhang }
1394814e85d6SHong Zhang 
139514d0ab18SJacob Faibussowitsch /*@C
139614d0ab18SJacob Faibussowitsch   TSAdjointSetFromOptions - Sets various `TS` adjoint parameters from options database.
1397814e85d6SHong Zhang 
1398c3339decSBarry Smith   Collective
1399814e85d6SHong Zhang 
140014d0ab18SJacob Faibussowitsch   Input Parameters:
140114d0ab18SJacob Faibussowitsch + ts                 - the `TS` context
140214d0ab18SJacob Faibussowitsch - PetscOptionsObject - the options context
1403814e85d6SHong Zhang 
1404814e85d6SHong Zhang   Options Database Keys:
140514d0ab18SJacob Faibussowitsch + -ts_adjoint_solve <yes,no>     - After solving the ODE/DAE solve the adjoint problem (requires `-ts_save_trajectory`)
1406814e85d6SHong Zhang . -ts_adjoint_monitor            - print information at each adjoint time step
1407814e85d6SHong Zhang - -ts_adjoint_monitor_draw_sensi - monitor the sensitivity of the first cost function wrt initial conditions (lambda[0]) graphically
1408814e85d6SHong Zhang 
1409814e85d6SHong Zhang   Level: developer
1410814e85d6SHong Zhang 
1411bcf0153eSBarry Smith   Note:
1412814e85d6SHong Zhang   This is not normally called directly by users
1413814e85d6SHong Zhang 
14141cc06b55SBarry Smith .seealso: [](ch_ts), `TSSetSaveTrajectory()`, `TSTrajectorySetUp()`
141514d0ab18SJacob Faibussowitsch @*/
1416d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointSetFromOptions(TS ts, PetscOptionItems *PetscOptionsObject)
1417d71ae5a4SJacob Faibussowitsch {
1418814e85d6SHong Zhang   PetscBool tflg, opt;
1419814e85d6SHong Zhang 
1420814e85d6SHong Zhang   PetscFunctionBegin;
1421dbbe0bcdSBarry Smith   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1422d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "TS Adjoint options");
1423814e85d6SHong Zhang   tflg = ts->adjoint_solve ? PETSC_TRUE : PETSC_FALSE;
14249566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-ts_adjoint_solve", "Solve the adjoint problem immediately after solving the forward problem", "", tflg, &tflg, &opt));
1425814e85d6SHong Zhang   if (opt) {
14269566063dSJacob Faibussowitsch     PetscCall(TSSetSaveTrajectory(ts));
1427814e85d6SHong Zhang     ts->adjoint_solve = tflg;
1428814e85d6SHong Zhang   }
14299566063dSJacob Faibussowitsch   PetscCall(TSAdjointMonitorSetFromOptions(ts, "-ts_adjoint_monitor", "Monitor adjoint timestep size", "TSAdjointMonitorDefault", TSAdjointMonitorDefault, NULL));
14309566063dSJacob Faibussowitsch   PetscCall(TSAdjointMonitorSetFromOptions(ts, "-ts_adjoint_monitor_sensi", "Monitor sensitivity in the adjoint computation", "TSAdjointMonitorSensi", TSAdjointMonitorSensi, NULL));
1431814e85d6SHong Zhang   opt = PETSC_FALSE;
14329566063dSJacob Faibussowitsch   PetscCall(PetscOptionsName("-ts_adjoint_monitor_draw_sensi", "Monitor adjoint sensitivities (lambda only) graphically", "TSAdjointMonitorDrawSensi", &opt));
1433814e85d6SHong Zhang   if (opt) {
1434814e85d6SHong Zhang     TSMonitorDrawCtx ctx;
1435814e85d6SHong Zhang     PetscInt         howoften = 1;
1436814e85d6SHong Zhang 
14379566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-ts_adjoint_monitor_draw_sensi", "Monitor adjoint sensitivities (lambda only) graphically", "TSAdjointMonitorDrawSensi", howoften, &howoften, NULL));
14389566063dSJacob Faibussowitsch     PetscCall(TSMonitorDrawCtxCreate(PetscObjectComm((PetscObject)ts), NULL, NULL, PETSC_DECIDE, PETSC_DECIDE, 300, 300, howoften, &ctx));
14399566063dSJacob Faibussowitsch     PetscCall(TSAdjointMonitorSet(ts, TSAdjointMonitorDrawSensi, ctx, (PetscErrorCode(*)(void **))TSMonitorDrawCtxDestroy));
1440814e85d6SHong Zhang   }
14413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1442814e85d6SHong Zhang }
1443814e85d6SHong Zhang 
1444814e85d6SHong Zhang /*@
1445814e85d6SHong Zhang   TSAdjointStep - Steps one time step backward in the adjoint run
1446814e85d6SHong Zhang 
1447c3339decSBarry Smith   Collective
1448814e85d6SHong Zhang 
1449814e85d6SHong Zhang   Input Parameter:
1450bcf0153eSBarry Smith . ts - the `TS` context obtained from `TSCreate()`
1451814e85d6SHong Zhang 
1452814e85d6SHong Zhang   Level: intermediate
1453814e85d6SHong Zhang 
14541cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdjointSetUp()`, `TSAdjointSolve()`
1455814e85d6SHong Zhang @*/
1456d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointStep(TS ts)
1457d71ae5a4SJacob Faibussowitsch {
1458814e85d6SHong Zhang   DM dm;
1459814e85d6SHong Zhang 
1460814e85d6SHong Zhang   PetscFunctionBegin;
1461814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
14629566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
14639566063dSJacob Faibussowitsch   PetscCall(TSAdjointSetUp(ts));
14647b0e2f17SHong Zhang   ts->steps--; /* must decrease the step index before the adjoint step is taken. */
1465814e85d6SHong Zhang 
1466814e85d6SHong Zhang   ts->reason     = TS_CONVERGED_ITERATING;
1467814e85d6SHong Zhang   ts->ptime_prev = ts->ptime;
14689566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(TS_AdjointStep, ts, 0, 0, 0));
1469dbbe0bcdSBarry Smith   PetscUseTypeMethod(ts, adjointstep);
14709566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(TS_AdjointStep, ts, 0, 0, 0));
14717b0e2f17SHong Zhang   ts->adjoint_steps++;
1472814e85d6SHong Zhang 
1473814e85d6SHong Zhang   if (ts->reason < 0) {
14743c633725SBarry Smith     PetscCheck(!ts->errorifstepfailed, PetscObjectComm((PetscObject)ts), PETSC_ERR_NOT_CONVERGED, "TSAdjointStep has failed due to %s", TSConvergedReasons[ts->reason]);
1475814e85d6SHong Zhang   } else if (!ts->reason) {
1476814e85d6SHong Zhang     if (ts->adjoint_steps >= ts->adjoint_max_steps) ts->reason = TS_CONVERGED_ITS;
1477814e85d6SHong Zhang   }
14783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1479814e85d6SHong Zhang }
1480814e85d6SHong Zhang 
1481814e85d6SHong Zhang /*@
1482814e85d6SHong Zhang   TSAdjointSolve - Solves the discrete ajoint problem for an ODE/DAE
1483814e85d6SHong Zhang 
1484c3339decSBarry Smith   Collective
1485bcf0153eSBarry Smith   `
1486b43aa488SJacob Faibussowitsch 
1487814e85d6SHong Zhang   Input Parameter:
1488bcf0153eSBarry Smith . ts - the `TS` context obtained from `TSCreate()`
1489814e85d6SHong Zhang 
1490bcf0153eSBarry Smith   Options Database Key:
1491814e85d6SHong Zhang . -ts_adjoint_view_solution <viewerinfo> - views the first gradient with respect to the initial values
1492814e85d6SHong Zhang 
1493814e85d6SHong Zhang   Level: intermediate
1494814e85d6SHong Zhang 
1495814e85d6SHong Zhang   Notes:
1496bcf0153eSBarry Smith   This must be called after a call to `TSSolve()` that solves the forward problem
1497814e85d6SHong Zhang 
1498bcf0153eSBarry Smith   By default this will integrate back to the initial time, one can use `TSAdjointSetSteps()` to step back to a later time
1499814e85d6SHong Zhang 
150042747ad1SJacob Faibussowitsch .seealso: [](ch_ts), `TSCreate()`, `TSSetCostGradients()`, `TSSetSolution()`, `TSAdjointStep()`
1501814e85d6SHong Zhang @*/
1502d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointSolve(TS ts)
1503d71ae5a4SJacob Faibussowitsch {
1504f1d62c27SHong Zhang   static PetscBool cite = PETSC_FALSE;
15057f59fb53SHong Zhang #if defined(TSADJOINT_STAGE)
15067f59fb53SHong Zhang   PetscLogStage adjoint_stage;
15077f59fb53SHong Zhang #endif
1508814e85d6SHong Zhang 
1509814e85d6SHong Zhang   PetscFunctionBegin;
1510814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1511421b783aSHong Zhang   PetscCall(PetscCitationsRegister("@article{Zhang2022tsadjoint,\n"
1512421b783aSHong Zhang                                    "  title         = {{PETSc TSAdjoint: A Discrete Adjoint ODE Solver for First-Order and Second-Order Sensitivity Analysis}},\n"
1513f1d62c27SHong Zhang                                    "  author        = {Zhang, Hong and Constantinescu, Emil M.  and Smith, Barry F.},\n"
1514421b783aSHong Zhang                                    "  journal       = {SIAM Journal on Scientific Computing},\n"
1515421b783aSHong Zhang                                    "  volume        = {44},\n"
1516421b783aSHong Zhang                                    "  number        = {1},\n"
1517421b783aSHong Zhang                                    "  pages         = {C1-C24},\n"
1518421b783aSHong Zhang                                    "  doi           = {10.1137/21M140078X},\n"
15199371c9d4SSatish Balay                                    "  year          = {2022}\n}\n",
15209371c9d4SSatish Balay                                    &cite));
15217f59fb53SHong Zhang #if defined(TSADJOINT_STAGE)
15229566063dSJacob Faibussowitsch   PetscCall(PetscLogStageRegister("TSAdjoint", &adjoint_stage));
15239566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePush(adjoint_stage));
15247f59fb53SHong Zhang #endif
15259566063dSJacob Faibussowitsch   PetscCall(TSAdjointSetUp(ts));
1526814e85d6SHong Zhang 
1527814e85d6SHong Zhang   /* reset time step and iteration counters */
1528814e85d6SHong Zhang   ts->adjoint_steps     = 0;
1529814e85d6SHong Zhang   ts->ksp_its           = 0;
1530814e85d6SHong Zhang   ts->snes_its          = 0;
1531814e85d6SHong Zhang   ts->num_snes_failures = 0;
1532814e85d6SHong Zhang   ts->reject            = 0;
1533814e85d6SHong Zhang   ts->reason            = TS_CONVERGED_ITERATING;
1534814e85d6SHong Zhang 
1535814e85d6SHong Zhang   if (!ts->adjoint_max_steps) ts->adjoint_max_steps = ts->steps;
1536814e85d6SHong Zhang   if (ts->adjoint_steps >= ts->adjoint_max_steps) ts->reason = TS_CONVERGED_ITS;
1537814e85d6SHong Zhang 
1538814e85d6SHong Zhang   while (!ts->reason) {
15399566063dSJacob Faibussowitsch     PetscCall(TSTrajectoryGet(ts->trajectory, ts, ts->steps, &ts->ptime));
15409566063dSJacob Faibussowitsch     PetscCall(TSAdjointMonitor(ts, ts->steps, ts->ptime, ts->vec_sol, ts->numcost, ts->vecs_sensi, ts->vecs_sensip));
15419566063dSJacob Faibussowitsch     PetscCall(TSAdjointEventHandler(ts));
15429566063dSJacob Faibussowitsch     PetscCall(TSAdjointStep(ts));
154348a46eb9SPierre Jolivet     if ((ts->vec_costintegral || ts->quadraturets) && !ts->costintegralfwd) PetscCall(TSAdjointCostIntegral(ts));
1544814e85d6SHong Zhang   }
1545bdbff044SHong Zhang   if (!ts->steps) {
15469566063dSJacob Faibussowitsch     PetscCall(TSTrajectoryGet(ts->trajectory, ts, ts->steps, &ts->ptime));
15479566063dSJacob Faibussowitsch     PetscCall(TSAdjointMonitor(ts, ts->steps, ts->ptime, ts->vec_sol, ts->numcost, ts->vecs_sensi, ts->vecs_sensip));
1548bdbff044SHong Zhang   }
1549814e85d6SHong Zhang   ts->solvetime = ts->ptime;
15509566063dSJacob Faibussowitsch   PetscCall(TSTrajectoryViewFromOptions(ts->trajectory, NULL, "-ts_trajectory_view"));
15519566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(ts->vecs_sensi[0], (PetscObject)ts, "-ts_adjoint_view_solution"));
1552814e85d6SHong Zhang   ts->adjoint_max_steps = 0;
15537f59fb53SHong Zhang #if defined(TSADJOINT_STAGE)
15549566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePop());
15557f59fb53SHong Zhang #endif
15563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1557814e85d6SHong Zhang }
1558814e85d6SHong Zhang 
1559814e85d6SHong Zhang /*@C
1560bcf0153eSBarry Smith   TSAdjointMonitor - Runs all user-provided adjoint monitor routines set using `TSAdjointMonitorSet()`
1561814e85d6SHong Zhang 
1562c3339decSBarry Smith   Collective
1563814e85d6SHong Zhang 
1564814e85d6SHong Zhang   Input Parameters:
1565bcf0153eSBarry Smith + ts      - time stepping context obtained from `TSCreate()`
1566814e85d6SHong Zhang . step    - step number that has just completed
1567814e85d6SHong Zhang . ptime   - model time of the state
1568814e85d6SHong Zhang . u       - state at the current model time
1569814e85d6SHong Zhang . numcost - number of cost functions (dimension of lambda  or mu)
1570814e85d6SHong Zhang . lambda  - vectors containing the gradients of the cost functions with respect to the ODE/DAE solution variables
1571814e85d6SHong Zhang - mu      - vectors containing the gradients of the cost functions with respect to the problem parameters
1572814e85d6SHong Zhang 
1573814e85d6SHong Zhang   Level: developer
1574814e85d6SHong Zhang 
1575bcf0153eSBarry Smith   Note:
1576bcf0153eSBarry Smith   `TSAdjointMonitor()` is typically used automatically within the time stepping implementations.
1577bcf0153eSBarry Smith   Users would almost never call this routine directly.
1578bcf0153eSBarry Smith 
1579bcf0153eSBarry Smith .seealso: `TSAdjointMonitorSet()`, `TSAdjointSolve()`
1580814e85d6SHong Zhang @*/
1581d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointMonitor(TS ts, PetscInt step, PetscReal ptime, Vec u, PetscInt numcost, Vec *lambda, Vec *mu)
1582d71ae5a4SJacob Faibussowitsch {
1583814e85d6SHong Zhang   PetscInt i, n = ts->numberadjointmonitors;
1584814e85d6SHong Zhang 
1585814e85d6SHong Zhang   PetscFunctionBegin;
1586814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1587814e85d6SHong Zhang   PetscValidHeaderSpecific(u, VEC_CLASSID, 4);
15889566063dSJacob Faibussowitsch   PetscCall(VecLockReadPush(u));
158948a46eb9SPierre Jolivet   for (i = 0; i < n; i++) PetscCall((*ts->adjointmonitor[i])(ts, step, ptime, u, numcost, lambda, mu, ts->adjointmonitorcontext[i]));
15909566063dSJacob Faibussowitsch   PetscCall(VecLockReadPop(u));
15913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1592814e85d6SHong Zhang }
1593814e85d6SHong Zhang 
1594814e85d6SHong Zhang /*@
1595814e85d6SHong Zhang   TSAdjointCostIntegral - Evaluate the cost integral in the adjoint run.
1596814e85d6SHong Zhang 
1597c3339decSBarry Smith   Collective
1598814e85d6SHong Zhang 
15994165533cSJose E. Roman   Input Parameter:
1600814e85d6SHong Zhang . ts - time stepping context
1601814e85d6SHong Zhang 
1602814e85d6SHong Zhang   Level: advanced
1603814e85d6SHong Zhang 
1604814e85d6SHong Zhang   Notes:
1605bcf0153eSBarry Smith   This function cannot be called until `TSAdjointStep()` has been completed.
1606814e85d6SHong Zhang 
16071cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdjointSolve()`, `TSAdjointStep()`
1608814e85d6SHong Zhang  @*/
1609d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointCostIntegral(TS ts)
1610d71ae5a4SJacob Faibussowitsch {
1611362febeeSStefano Zampini   PetscFunctionBegin;
1612814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1613dbbe0bcdSBarry Smith   PetscUseTypeMethod(ts, adjointintegral);
16143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1615814e85d6SHong Zhang }
1616814e85d6SHong Zhang 
1617814e85d6SHong Zhang /* ------------------ Forward (tangent linear) sensitivity  ------------------*/
1618814e85d6SHong Zhang 
1619814e85d6SHong Zhang /*@
1620814e85d6SHong Zhang   TSForwardSetUp - Sets up the internal data structures for the later use
1621814e85d6SHong Zhang   of forward sensitivity analysis
1622814e85d6SHong Zhang 
1623c3339decSBarry Smith   Collective
1624814e85d6SHong Zhang 
1625814e85d6SHong Zhang   Input Parameter:
1626bcf0153eSBarry Smith . ts - the `TS` context obtained from `TSCreate()`
1627814e85d6SHong Zhang 
1628814e85d6SHong Zhang   Level: advanced
1629814e85d6SHong Zhang 
16301cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSCreate()`, `TSDestroy()`, `TSSetUp()`
1631814e85d6SHong Zhang @*/
1632d71ae5a4SJacob Faibussowitsch PetscErrorCode TSForwardSetUp(TS ts)
1633d71ae5a4SJacob Faibussowitsch {
1634814e85d6SHong Zhang   PetscFunctionBegin;
1635814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
16363ba16761SJacob Faibussowitsch   if (ts->forwardsetupcalled) PetscFunctionReturn(PETSC_SUCCESS);
1637dbbe0bcdSBarry Smith   PetscTryTypeMethod(ts, forwardsetup);
16389566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &ts->vec_sensip_col));
1639814e85d6SHong Zhang   ts->forwardsetupcalled = PETSC_TRUE;
16403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1641814e85d6SHong Zhang }
1642814e85d6SHong Zhang 
1643814e85d6SHong Zhang /*@
16447adebddeSHong Zhang   TSForwardReset - Reset the internal data structures used by forward sensitivity analysis
16457adebddeSHong Zhang 
1646c3339decSBarry Smith   Collective
16477adebddeSHong Zhang 
16487adebddeSHong Zhang   Input Parameter:
1649bcf0153eSBarry Smith . ts - the `TS` context obtained from `TSCreate()`
16507adebddeSHong Zhang 
16517adebddeSHong Zhang   Level: advanced
16527adebddeSHong Zhang 
16531cc06b55SBarry Smith .seealso: [](ch_ts), `TSCreate()`, `TSDestroy()`, `TSForwardSetUp()`
16547adebddeSHong Zhang @*/
1655d71ae5a4SJacob Faibussowitsch PetscErrorCode TSForwardReset(TS ts)
1656d71ae5a4SJacob Faibussowitsch {
16577207e0fdSHong Zhang   TS quadts = ts->quadraturets;
16587adebddeSHong Zhang 
16597adebddeSHong Zhang   PetscFunctionBegin;
16607adebddeSHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1661dbbe0bcdSBarry Smith   PetscTryTypeMethod(ts, forwardreset);
16629566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&ts->mat_sensip));
166348a46eb9SPierre Jolivet   if (quadts) PetscCall(MatDestroy(&quadts->mat_sensip));
16649566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ts->vec_sensip_col));
16657adebddeSHong Zhang   ts->forward_solve      = PETSC_FALSE;
16667adebddeSHong Zhang   ts->forwardsetupcalled = PETSC_FALSE;
16673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16687adebddeSHong Zhang }
16697adebddeSHong Zhang 
16707adebddeSHong Zhang /*@
1671814e85d6SHong Zhang   TSForwardSetIntegralGradients - Set the vectors holding forward sensitivities of the integral term.
1672814e85d6SHong Zhang 
1673d8d19677SJose E. Roman   Input Parameters:
1674bcf0153eSBarry Smith + ts        - the `TS` context obtained from `TSCreate()`
1675814e85d6SHong Zhang . numfwdint - number of integrals
167667b8a455SSatish Balay - vp        - the vectors containing the gradients for each integral w.r.t. parameters
1677814e85d6SHong Zhang 
16787207e0fdSHong Zhang   Level: deprecated
1679814e85d6SHong Zhang 
168042747ad1SJacob Faibussowitsch .seealso: [](ch_ts), `TSForwardGetSensitivities()`, `TSForwardGetIntegralGradients()`, `TSForwardStep()`
1681814e85d6SHong Zhang @*/
1682d71ae5a4SJacob Faibussowitsch PetscErrorCode TSForwardSetIntegralGradients(TS ts, PetscInt numfwdint, Vec *vp)
1683d71ae5a4SJacob Faibussowitsch {
1684814e85d6SHong Zhang   PetscFunctionBegin;
1685814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
16863c633725SBarry 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()");
1687814e85d6SHong Zhang   if (!ts->numcost) ts->numcost = numfwdint;
1688814e85d6SHong Zhang 
1689814e85d6SHong Zhang   ts->vecs_integral_sensip = vp;
16903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1691814e85d6SHong Zhang }
1692814e85d6SHong Zhang 
1693814e85d6SHong Zhang /*@
1694814e85d6SHong Zhang   TSForwardGetIntegralGradients - Returns the forward sensitivities of the integral term.
1695814e85d6SHong Zhang 
1696814e85d6SHong Zhang   Input Parameter:
1697bcf0153eSBarry Smith . ts - the `TS` context obtained from `TSCreate()`
1698814e85d6SHong Zhang 
169914d0ab18SJacob Faibussowitsch   Output Parameters:
170014d0ab18SJacob Faibussowitsch + numfwdint - number of integrals
170114d0ab18SJacob Faibussowitsch - vp        - the vectors containing the gradients for each integral w.r.t. parameters
1702814e85d6SHong Zhang 
17037207e0fdSHong Zhang   Level: deprecated
1704814e85d6SHong Zhang 
170542747ad1SJacob Faibussowitsch .seealso: [](ch_ts), `TSForwardSetSensitivities()`, `TSForwardSetIntegralGradients()`, `TSForwardStep()`
1706814e85d6SHong Zhang @*/
1707d71ae5a4SJacob Faibussowitsch PetscErrorCode TSForwardGetIntegralGradients(TS ts, PetscInt *numfwdint, Vec **vp)
1708d71ae5a4SJacob Faibussowitsch {
1709814e85d6SHong Zhang   PetscFunctionBegin;
1710814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
17114f572ea9SToby Isaac   PetscAssertPointer(vp, 3);
1712814e85d6SHong Zhang   if (numfwdint) *numfwdint = ts->numcost;
1713814e85d6SHong Zhang   if (vp) *vp = ts->vecs_integral_sensip;
17143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1715814e85d6SHong Zhang }
1716814e85d6SHong Zhang 
1717814e85d6SHong Zhang /*@
1718814e85d6SHong Zhang   TSForwardStep - Compute the forward sensitivity for one time step.
1719814e85d6SHong Zhang 
1720c3339decSBarry Smith   Collective
1721814e85d6SHong Zhang 
17224165533cSJose E. Roman   Input Parameter:
1723814e85d6SHong Zhang . ts - time stepping context
1724814e85d6SHong Zhang 
1725814e85d6SHong Zhang   Level: advanced
1726814e85d6SHong Zhang 
1727814e85d6SHong Zhang   Notes:
1728bcf0153eSBarry Smith   This function cannot be called until `TSStep()` has been completed.
1729814e85d6SHong Zhang 
17301cc06b55SBarry Smith .seealso: [](ch_ts), `TSForwardSetSensitivities()`, `TSForwardGetSensitivities()`, `TSForwardSetIntegralGradients()`, `TSForwardGetIntegralGradients()`, `TSForwardSetUp()`
1731814e85d6SHong Zhang @*/
1732d71ae5a4SJacob Faibussowitsch PetscErrorCode TSForwardStep(TS ts)
1733d71ae5a4SJacob Faibussowitsch {
1734362febeeSStefano Zampini   PetscFunctionBegin;
1735814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
17369566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(TS_ForwardStep, ts, 0, 0, 0));
1737dbbe0bcdSBarry Smith   PetscUseTypeMethod(ts, forwardstep);
17389566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(TS_ForwardStep, ts, 0, 0, 0));
17393c633725SBarry Smith   PetscCheck(ts->reason >= 0 || !ts->errorifstepfailed, PetscObjectComm((PetscObject)ts), PETSC_ERR_NOT_CONVERGED, "TSFowardStep has failed due to %s", TSConvergedReasons[ts->reason]);
17403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1741814e85d6SHong Zhang }
1742814e85d6SHong Zhang 
1743814e85d6SHong Zhang /*@
1744814e85d6SHong Zhang   TSForwardSetSensitivities - Sets the initial value of the trajectory sensitivities of solution  w.r.t. the problem parameters and initial values.
1745814e85d6SHong Zhang 
1746c3339decSBarry Smith   Logically Collective
1747814e85d6SHong Zhang 
1748814e85d6SHong Zhang   Input Parameters:
1749bcf0153eSBarry Smith + ts   - the `TS` context obtained from `TSCreate()`
1750814e85d6SHong Zhang . nump - number of parameters
1751814e85d6SHong Zhang - Smat - sensitivities with respect to the parameters, the number of entries in these vectors is the same as the number of parameters
1752814e85d6SHong Zhang 
1753814e85d6SHong Zhang   Level: beginner
1754814e85d6SHong Zhang 
1755814e85d6SHong Zhang   Notes:
1756814e85d6SHong Zhang   Forward sensitivity is also called 'trajectory sensitivity' in some fields such as power systems.
1757bcf0153eSBarry Smith   This function turns on a flag to trigger `TSSolve()` to compute forward sensitivities automatically.
1758bcf0153eSBarry Smith   You must call this function before `TSSolve()`.
1759814e85d6SHong Zhang   The entries in the sensitivity matrix must be correctly initialized with the values S = dy/dp|startingtime.
1760814e85d6SHong Zhang 
17611cc06b55SBarry Smith .seealso: [](ch_ts), `TSForwardGetSensitivities()`, `TSForwardSetIntegralGradients()`, `TSForwardGetIntegralGradients()`, `TSForwardStep()`
1762814e85d6SHong Zhang @*/
1763d71ae5a4SJacob Faibussowitsch PetscErrorCode TSForwardSetSensitivities(TS ts, PetscInt nump, Mat Smat)
1764d71ae5a4SJacob Faibussowitsch {
1765814e85d6SHong Zhang   PetscFunctionBegin;
1766814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1767814e85d6SHong Zhang   PetscValidHeaderSpecific(Smat, MAT_CLASSID, 3);
1768814e85d6SHong Zhang   ts->forward_solve = PETSC_TRUE;
1769b5b4017aSHong Zhang   if (nump == PETSC_DEFAULT) {
17709566063dSJacob Faibussowitsch     PetscCall(MatGetSize(Smat, NULL, &ts->num_parameters));
1771b5b4017aSHong Zhang   } else ts->num_parameters = nump;
17729566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)Smat));
17739566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&ts->mat_sensip));
1774814e85d6SHong Zhang   ts->mat_sensip = Smat;
17753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1776814e85d6SHong Zhang }
1777814e85d6SHong Zhang 
1778814e85d6SHong Zhang /*@
1779814e85d6SHong Zhang   TSForwardGetSensitivities - Returns the trajectory sensitivities
1780814e85d6SHong Zhang 
1781bcf0153eSBarry Smith   Not Collective, but Smat returned is parallel if ts is parallel
1782814e85d6SHong Zhang 
1783d8d19677SJose E. Roman   Output Parameters:
1784bcf0153eSBarry Smith + ts   - the `TS` context obtained from `TSCreate()`
1785814e85d6SHong Zhang . nump - number of parameters
1786814e85d6SHong Zhang - Smat - sensitivities with respect to the parameters, the number of entries in these vectors is the same as the number of parameters
1787814e85d6SHong Zhang 
1788814e85d6SHong Zhang   Level: intermediate
1789814e85d6SHong Zhang 
17901cc06b55SBarry Smith .seealso: [](ch_ts), `TSForwardSetSensitivities()`, `TSForwardSetIntegralGradients()`, `TSForwardGetIntegralGradients()`, `TSForwardStep()`
1791814e85d6SHong Zhang @*/
1792d71ae5a4SJacob Faibussowitsch PetscErrorCode TSForwardGetSensitivities(TS ts, PetscInt *nump, Mat *Smat)
1793d71ae5a4SJacob Faibussowitsch {
1794814e85d6SHong Zhang   PetscFunctionBegin;
1795814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1796814e85d6SHong Zhang   if (nump) *nump = ts->num_parameters;
1797814e85d6SHong Zhang   if (Smat) *Smat = ts->mat_sensip;
17983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1799814e85d6SHong Zhang }
1800814e85d6SHong Zhang 
1801814e85d6SHong Zhang /*@
1802814e85d6SHong Zhang   TSForwardCostIntegral - Evaluate the cost integral in the forward run.
1803814e85d6SHong Zhang 
1804c3339decSBarry Smith   Collective
1805814e85d6SHong Zhang 
18064165533cSJose E. Roman   Input Parameter:
1807814e85d6SHong Zhang . ts - time stepping context
1808814e85d6SHong Zhang 
1809814e85d6SHong Zhang   Level: advanced
1810814e85d6SHong Zhang 
1811bcf0153eSBarry Smith   Note:
1812bcf0153eSBarry Smith   This function cannot be called until `TSStep()` has been completed.
1813814e85d6SHong Zhang 
18141cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSAdjointCostIntegral()`
1815814e85d6SHong Zhang @*/
1816d71ae5a4SJacob Faibussowitsch PetscErrorCode TSForwardCostIntegral(TS ts)
1817d71ae5a4SJacob Faibussowitsch {
1818362febeeSStefano Zampini   PetscFunctionBegin;
1819814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1820dbbe0bcdSBarry Smith   PetscUseTypeMethod(ts, forwardintegral);
18213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1822814e85d6SHong Zhang }
1823b5b4017aSHong Zhang 
1824b5b4017aSHong Zhang /*@
1825b5b4017aSHong Zhang   TSForwardSetInitialSensitivities - Set initial values for tangent linear sensitivities
1826b5b4017aSHong Zhang 
1827c3339decSBarry Smith   Collective
1828b5b4017aSHong Zhang 
1829d8d19677SJose E. Roman   Input Parameters:
1830bcf0153eSBarry Smith + ts   - the `TS` context obtained from `TSCreate()`
1831b5b4017aSHong Zhang - didp - parametric sensitivities of the initial condition
1832b5b4017aSHong Zhang 
1833b5b4017aSHong Zhang   Level: intermediate
1834b5b4017aSHong Zhang 
1835bcf0153eSBarry Smith   Notes:
1836bcf0153eSBarry Smith   `TSSolve()` allows users to pass the initial solution directly to `TS`. But the tangent linear variables cannot be initialized in this way.
1837bcf0153eSBarry Smith   This function is used to set initial values for tangent linear variables.
1838b5b4017aSHong Zhang 
18391cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSForwardSetSensitivities()`
1840b5b4017aSHong Zhang @*/
1841d71ae5a4SJacob Faibussowitsch PetscErrorCode TSForwardSetInitialSensitivities(TS ts, Mat didp)
1842d71ae5a4SJacob Faibussowitsch {
1843362febeeSStefano Zampini   PetscFunctionBegin;
1844b5b4017aSHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1845b5b4017aSHong Zhang   PetscValidHeaderSpecific(didp, MAT_CLASSID, 2);
184648a46eb9SPierre Jolivet   if (!ts->mat_sensip) PetscCall(TSForwardSetSensitivities(ts, PETSC_DEFAULT, didp));
18473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1848b5b4017aSHong Zhang }
18496affb6f8SHong Zhang 
18506affb6f8SHong Zhang /*@
18516affb6f8SHong Zhang   TSForwardGetStages - Get the number of stages and the tangent linear sensitivities at the intermediate stages
18526affb6f8SHong Zhang 
18536affb6f8SHong Zhang   Input Parameter:
1854bcf0153eSBarry Smith . ts - the `TS` context obtained from `TSCreate()`
18556affb6f8SHong Zhang 
18566affb6f8SHong Zhang   Output Parameters:
1857cd4cee2dSHong Zhang + ns - number of stages
18586affb6f8SHong Zhang - S  - tangent linear sensitivities at the intermediate stages
18596affb6f8SHong Zhang 
18606affb6f8SHong Zhang   Level: advanced
18616affb6f8SHong Zhang 
1862bcf0153eSBarry Smith .seealso: `TS`
18636affb6f8SHong Zhang @*/
1864d71ae5a4SJacob Faibussowitsch PetscErrorCode TSForwardGetStages(TS ts, PetscInt *ns, Mat **S)
1865d71ae5a4SJacob Faibussowitsch {
18666affb6f8SHong Zhang   PetscFunctionBegin;
18676affb6f8SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
18686affb6f8SHong Zhang 
18696affb6f8SHong Zhang   if (!ts->ops->getstages) *S = NULL;
1870dbbe0bcdSBarry Smith   else PetscUseTypeMethod(ts, forwardgetstages, ns, S);
18713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
18726affb6f8SHong Zhang }
1873cd4cee2dSHong Zhang 
1874cd4cee2dSHong Zhang /*@
1875bcf0153eSBarry Smith   TSCreateQuadratureTS - Create a sub-`TS` that evaluates integrals over time
1876cd4cee2dSHong Zhang 
1877d8d19677SJose E. Roman   Input Parameters:
1878bcf0153eSBarry Smith + ts  - the `TS` context obtained from `TSCreate()`
1879cd4cee2dSHong Zhang - fwd - flag indicating whether to evaluate cost integral in the forward run or the adjoint run
1880cd4cee2dSHong Zhang 
18812fe279fdSBarry Smith   Output Parameter:
1882bcf0153eSBarry Smith . quadts - the child `TS` context
1883cd4cee2dSHong Zhang 
1884cd4cee2dSHong Zhang   Level: intermediate
1885cd4cee2dSHong Zhang 
18861cc06b55SBarry Smith .seealso: [](ch_ts), `TSGetQuadratureTS()`
1887cd4cee2dSHong Zhang @*/
1888d71ae5a4SJacob Faibussowitsch PetscErrorCode TSCreateQuadratureTS(TS ts, PetscBool fwd, TS *quadts)
1889d71ae5a4SJacob Faibussowitsch {
1890cd4cee2dSHong Zhang   char prefix[128];
1891cd4cee2dSHong Zhang 
1892cd4cee2dSHong Zhang   PetscFunctionBegin;
1893cd4cee2dSHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
18944f572ea9SToby Isaac   PetscAssertPointer(quadts, 3);
18959566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts->quadraturets));
18969566063dSJacob Faibussowitsch   PetscCall(TSCreate(PetscObjectComm((PetscObject)ts), &ts->quadraturets));
18979566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)ts->quadraturets, (PetscObject)ts, 1));
18989566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(prefix, sizeof(prefix), "%squad_", ((PetscObject)ts)->prefix ? ((PetscObject)ts)->prefix : ""));
18999566063dSJacob Faibussowitsch   PetscCall(TSSetOptionsPrefix(ts->quadraturets, prefix));
1900cd4cee2dSHong Zhang   *quadts = ts->quadraturets;
1901cd4cee2dSHong Zhang 
1902cd4cee2dSHong Zhang   if (ts->numcost) {
19039566063dSJacob Faibussowitsch     PetscCall(VecCreateSeq(PETSC_COMM_SELF, ts->numcost, &(*quadts)->vec_sol));
1904cd4cee2dSHong Zhang   } else {
19059566063dSJacob Faibussowitsch     PetscCall(VecCreateSeq(PETSC_COMM_SELF, 1, &(*quadts)->vec_sol));
1906cd4cee2dSHong Zhang   }
1907cd4cee2dSHong Zhang   ts->costintegralfwd = fwd;
19083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1909cd4cee2dSHong Zhang }
1910cd4cee2dSHong Zhang 
1911cd4cee2dSHong Zhang /*@
1912bcf0153eSBarry Smith   TSGetQuadratureTS - Return the sub-`TS` that evaluates integrals over time
1913cd4cee2dSHong Zhang 
1914cd4cee2dSHong Zhang   Input Parameter:
1915bcf0153eSBarry Smith . ts - the `TS` context obtained from `TSCreate()`
1916cd4cee2dSHong Zhang 
1917cd4cee2dSHong Zhang   Output Parameters:
1918cd4cee2dSHong Zhang + fwd    - flag indicating whether to evaluate cost integral in the forward run or the adjoint run
1919bcf0153eSBarry Smith - quadts - the child `TS` context
1920cd4cee2dSHong Zhang 
1921cd4cee2dSHong Zhang   Level: intermediate
1922cd4cee2dSHong Zhang 
19231cc06b55SBarry Smith .seealso: [](ch_ts), `TSCreateQuadratureTS()`
1924cd4cee2dSHong Zhang @*/
1925d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGetQuadratureTS(TS ts, PetscBool *fwd, TS *quadts)
1926d71ae5a4SJacob Faibussowitsch {
1927cd4cee2dSHong Zhang   PetscFunctionBegin;
1928cd4cee2dSHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1929cd4cee2dSHong Zhang   if (fwd) *fwd = ts->costintegralfwd;
1930cd4cee2dSHong Zhang   if (quadts) *quadts = ts->quadraturets;
19313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1932cd4cee2dSHong Zhang }
1933ba0675f6SHong Zhang 
1934ba0675f6SHong Zhang /*@
1935bcf0153eSBarry Smith   TSComputeSNESJacobian - Compute the Jacobian needed for the `SNESSolve()` in `TS`
1936bcf0153eSBarry Smith 
1937c3339decSBarry Smith   Collective
1938ba0675f6SHong Zhang 
1939ba0675f6SHong Zhang   Input Parameters:
1940bcf0153eSBarry Smith + ts - the `TS` context obtained from `TSCreate()`
1941ba0675f6SHong Zhang - x  - state vector
1942ba0675f6SHong Zhang 
1943ba0675f6SHong Zhang   Output Parameters:
1944ba0675f6SHong Zhang + J    - Jacobian matrix
1945ba0675f6SHong Zhang - Jpre - preconditioning matrix for J (may be same as J)
1946ba0675f6SHong Zhang 
1947ba0675f6SHong Zhang   Level: developer
1948ba0675f6SHong Zhang 
1949bcf0153eSBarry Smith   Note:
1950bcf0153eSBarry Smith   Uses finite differencing when `TS` Jacobian is not available.
1951bcf0153eSBarry Smith 
1952b43aa488SJacob Faibussowitsch .seealso: `SNES`, `TS`, `SNESSetJacobian()`, `TSSetRHSJacobian()`, `TSSetIJacobian()`
1953ba0675f6SHong Zhang @*/
1954d71ae5a4SJacob Faibussowitsch PetscErrorCode TSComputeSNESJacobian(TS ts, Vec x, Mat J, Mat Jpre)
1955d71ae5a4SJacob Faibussowitsch {
1956ba0675f6SHong Zhang   SNES snes                                          = ts->snes;
1957ba0675f6SHong Zhang   PetscErrorCode (*jac)(SNES, Vec, Mat, Mat, void *) = NULL;
1958ba0675f6SHong Zhang 
1959ba0675f6SHong Zhang   PetscFunctionBegin;
1960ba0675f6SHong Zhang   /*
1961ba0675f6SHong Zhang     Unlike implicit methods, explicit methods do not have SNESMatFDColoring in the snes object
1962ba0675f6SHong Zhang     because SNESSolve() has not been called yet; so querying SNESMatFDColoring does not work for
1963da81f932SPierre Jolivet     explicit methods. Instead, we check the Jacobian compute function directly to determine if FD
1964ba0675f6SHong Zhang     coloring is used.
1965ba0675f6SHong Zhang   */
19669566063dSJacob Faibussowitsch   PetscCall(SNESGetJacobian(snes, NULL, NULL, &jac, NULL));
1967ba0675f6SHong Zhang   if (jac == SNESComputeJacobianDefaultColor) {
1968ba0675f6SHong Zhang     Vec f;
19699566063dSJacob Faibussowitsch     PetscCall(SNESSetSolution(snes, x));
19709566063dSJacob Faibussowitsch     PetscCall(SNESGetFunction(snes, &f, NULL, NULL));
1971ba0675f6SHong Zhang     /* Force MatFDColoringApply to evaluate the SNES residual function for the base vector */
19729566063dSJacob Faibussowitsch     PetscCall(SNESComputeFunction(snes, x, f));
1973ba0675f6SHong Zhang   }
19749566063dSJacob Faibussowitsch   PetscCall(SNESComputeJacobian(snes, x, J, Jpre));
19753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1976ba0675f6SHong Zhang }
1977