xref: /petsc/src/ts/interface/sensitivity/tssen.c (revision 33f72c5d48e626d09b7d7a978ff73ecaa099928a)
1814e85d6SHong Zhang #include <petsc/private/tsimpl.h>        /*I "petscts.h"  I*/
2814e85d6SHong Zhang #include <petscdraw.h>
3814e85d6SHong Zhang 
4*33f72c5dSHong Zhang PetscLogEvent TS_AdjointStep,TS_ForwardStep,TS_JacobianPEval;
5814e85d6SHong Zhang 
6814e85d6SHong Zhang /* ------------------------ Sensitivity Context ---------------------------*/
7814e85d6SHong Zhang 
8814e85d6SHong Zhang /*@C
9b5b4017aSHong Zhang   TSSetRHSJacobianP - Sets the function that computes the Jacobian of G w.r.t. the parameters P where U_t = G(U,P,t), as well as the location to store the matrix.
10814e85d6SHong Zhang 
11814e85d6SHong Zhang   Logically Collective on TS
12814e85d6SHong Zhang 
13814e85d6SHong Zhang   Input Parameters:
14b5b4017aSHong Zhang + ts - TS context obtained from TSCreate()
15b5b4017aSHong Zhang . Amat - JacobianP matrix
16b5b4017aSHong Zhang . func - function
17b5b4017aSHong Zhang - ctx - [optional] user-defined function context
18814e85d6SHong Zhang 
19814e85d6SHong Zhang   Calling sequence of func:
20814e85d6SHong Zhang $ func (TS ts,PetscReal t,Vec y,Mat A,void *ctx);
21814e85d6SHong Zhang +   t - current timestep
22c9ad9de0SHong Zhang .   U - input vector (current ODE solution)
23814e85d6SHong Zhang .   A - output matrix
24814e85d6SHong Zhang -   ctx - [optional] user-defined function context
25814e85d6SHong Zhang 
26814e85d6SHong Zhang   Level: intermediate
27814e85d6SHong Zhang 
28814e85d6SHong Zhang   Notes:
29814e85d6SHong Zhang     Amat has the same number of rows and the same row parallel layout as u, Amat has the same number of columns and parallel layout as p
30814e85d6SHong Zhang 
31814e85d6SHong Zhang .keywords: TS, sensitivity
32814e85d6SHong Zhang .seealso:
33814e85d6SHong Zhang @*/
34814e85d6SHong Zhang PetscErrorCode TSSetRHSJacobianP(TS ts,Mat Amat,PetscErrorCode (*func)(TS,PetscReal,Vec,Mat,void*),void *ctx)
35814e85d6SHong Zhang {
36814e85d6SHong Zhang   PetscErrorCode ierr;
37814e85d6SHong Zhang 
38814e85d6SHong Zhang   PetscFunctionBegin;
39814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
40814e85d6SHong Zhang   PetscValidHeaderSpecific(Amat,MAT_CLASSID,2);
41814e85d6SHong Zhang 
42814e85d6SHong Zhang   ts->rhsjacobianp    = func;
43814e85d6SHong Zhang   ts->rhsjacobianpctx = ctx;
44814e85d6SHong Zhang   if(Amat) {
45814e85d6SHong Zhang     ierr = PetscObjectReference((PetscObject)Amat);CHKERRQ(ierr);
46*33f72c5dSHong Zhang     ierr = MatDestroy(&ts->Jacprhs);CHKERRQ(ierr);
47*33f72c5dSHong Zhang     ts->Jacprhs = Amat;
48814e85d6SHong Zhang   }
49814e85d6SHong Zhang   PetscFunctionReturn(0);
50814e85d6SHong Zhang }
51814e85d6SHong Zhang 
52814e85d6SHong Zhang /*@C
53814e85d6SHong Zhang   TSComputeRHSJacobianP - Runs the user-defined JacobianP function.
54814e85d6SHong Zhang 
55814e85d6SHong Zhang   Collective on TS
56814e85d6SHong Zhang 
57814e85d6SHong Zhang   Input Parameters:
58814e85d6SHong Zhang . ts   - The TS context obtained from TSCreate()
59814e85d6SHong Zhang 
60814e85d6SHong Zhang   Level: developer
61814e85d6SHong Zhang 
62814e85d6SHong Zhang .keywords: TS, sensitivity
63814e85d6SHong Zhang .seealso: TSSetRHSJacobianP()
64814e85d6SHong Zhang @*/
65c9ad9de0SHong Zhang PetscErrorCode TSComputeRHSJacobianP(TS ts,PetscReal t,Vec U,Mat Amat)
66814e85d6SHong Zhang {
67814e85d6SHong Zhang   PetscErrorCode ierr;
68814e85d6SHong Zhang 
69814e85d6SHong Zhang   PetscFunctionBegin;
70*33f72c5dSHong Zhang   if (!Amat) PetscFunctionReturn(0);
71814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
72c9ad9de0SHong Zhang   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
73814e85d6SHong Zhang 
74814e85d6SHong Zhang   PetscStackPush("TS user JacobianP function for sensitivity analysis");
75c9ad9de0SHong Zhang   ierr = (*ts->rhsjacobianp)(ts,t,U,Amat,ts->rhsjacobianpctx);CHKERRQ(ierr);
76814e85d6SHong Zhang   PetscStackPop;
77814e85d6SHong Zhang   PetscFunctionReturn(0);
78814e85d6SHong Zhang }
79814e85d6SHong Zhang 
80814e85d6SHong Zhang /*@C
81*33f72c5dSHong 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.
82*33f72c5dSHong Zhang 
83*33f72c5dSHong Zhang   Logically Collective on TS
84*33f72c5dSHong Zhang 
85*33f72c5dSHong Zhang   Input Parameters:
86*33f72c5dSHong Zhang + ts - TS context obtained from TSCreate()
87*33f72c5dSHong Zhang . Amat - JacobianP matrix
88*33f72c5dSHong Zhang . func - function
89*33f72c5dSHong Zhang - ctx - [optional] user-defined function context
90*33f72c5dSHong Zhang 
91*33f72c5dSHong Zhang   Calling sequence of func:
92*33f72c5dSHong Zhang $ func (TS ts,PetscReal t,Vec y,Mat A,void *ctx);
93*33f72c5dSHong Zhang +   t - current timestep
94*33f72c5dSHong Zhang .   U - input vector (current ODE solution)
95*33f72c5dSHong Zhang .   Udot - time derivative of state vector
96*33f72c5dSHong Zhang .   shift - shift to apply, see note below
97*33f72c5dSHong Zhang .   A - output matrix
98*33f72c5dSHong Zhang -   ctx - [optional] user-defined function context
99*33f72c5dSHong Zhang 
100*33f72c5dSHong Zhang   Level: intermediate
101*33f72c5dSHong Zhang 
102*33f72c5dSHong Zhang   Notes:
103*33f72c5dSHong 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
104*33f72c5dSHong Zhang 
105*33f72c5dSHong Zhang .keywords: TS, sensitivity
106*33f72c5dSHong Zhang .seealso:
107*33f72c5dSHong Zhang @*/
108*33f72c5dSHong Zhang PetscErrorCode TSSetIJacobianP(TS ts,Mat Amat,PetscErrorCode (*func)(TS,PetscReal,Vec,Vec,PetscReal,Mat,void*),void *ctx)
109*33f72c5dSHong Zhang {
110*33f72c5dSHong Zhang   PetscErrorCode ierr;
111*33f72c5dSHong Zhang 
112*33f72c5dSHong Zhang   PetscFunctionBegin;
113*33f72c5dSHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
114*33f72c5dSHong Zhang   PetscValidHeaderSpecific(Amat,MAT_CLASSID,2);
115*33f72c5dSHong Zhang 
116*33f72c5dSHong Zhang   ts->ijacobianp    = func;
117*33f72c5dSHong Zhang   ts->ijacobianpctx = ctx;
118*33f72c5dSHong Zhang   if(Amat) {
119*33f72c5dSHong Zhang     ierr = PetscObjectReference((PetscObject)Amat);CHKERRQ(ierr);
120*33f72c5dSHong Zhang     ierr = MatDestroy(&ts->Jacp);CHKERRQ(ierr);
121*33f72c5dSHong Zhang     ts->Jacp = Amat;
122*33f72c5dSHong Zhang   }
123*33f72c5dSHong Zhang   PetscFunctionReturn(0);
124*33f72c5dSHong Zhang }
125*33f72c5dSHong Zhang 
126*33f72c5dSHong Zhang /*@C
127*33f72c5dSHong Zhang   TSComputeIJacobianP - Runs the user-defined IJacobianP function.
128*33f72c5dSHong Zhang 
129*33f72c5dSHong Zhang   Collective on TS
130*33f72c5dSHong Zhang 
131*33f72c5dSHong Zhang   Input Parameters:
132*33f72c5dSHong Zhang + ts - the TS context
133*33f72c5dSHong Zhang . t - current timestep
134*33f72c5dSHong Zhang . U - state vector
135*33f72c5dSHong Zhang . Udot - time derivative of state vector
136*33f72c5dSHong Zhang . shift - shift to apply, see note below
137*33f72c5dSHong Zhang - imex - flag indicates if the method is IMEX so that the RHSJacobian should be kept separate
138*33f72c5dSHong Zhang 
139*33f72c5dSHong Zhang   Output Parameters:
140*33f72c5dSHong Zhang . A - Jacobian matrix
141*33f72c5dSHong Zhang 
142*33f72c5dSHong Zhang   Level: developer
143*33f72c5dSHong Zhang 
144*33f72c5dSHong Zhang .keywords: TS, sensitivity
145*33f72c5dSHong Zhang .seealso: TSSetIJacobianP()
146*33f72c5dSHong Zhang @*/
147*33f72c5dSHong Zhang PetscErrorCode TSComputeIJacobianP(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal shift,Mat Amat,PetscBool imex)
148*33f72c5dSHong Zhang {
149*33f72c5dSHong Zhang   PetscErrorCode ierr;
150*33f72c5dSHong Zhang 
151*33f72c5dSHong Zhang   PetscFunctionBegin;
152*33f72c5dSHong Zhang   if (!Amat) PetscFunctionReturn(0);
153*33f72c5dSHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
154*33f72c5dSHong Zhang   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
155*33f72c5dSHong Zhang   PetscValidHeaderSpecific(Udot,VEC_CLASSID,4);
156*33f72c5dSHong Zhang 
157*33f72c5dSHong Zhang   ierr = PetscLogEventBegin(TS_JacobianPEval,ts,U,Amat,0);CHKERRQ(ierr);
158*33f72c5dSHong Zhang   if (ts->ijacobianp) {
159*33f72c5dSHong Zhang     PetscStackPush("TS user JacobianP function for sensitivity analysis");
160*33f72c5dSHong Zhang     ierr = (*ts->ijacobianp)(ts,t,U,Udot,shift,Amat,ts->ijacobianpctx);CHKERRQ(ierr);
161*33f72c5dSHong Zhang     PetscStackPop;
162*33f72c5dSHong Zhang   }
163*33f72c5dSHong Zhang   if (imex) {
164*33f72c5dSHong Zhang     if (!ts->ijacobianp) {  /* system was written as Udot = G(t,U) */
165*33f72c5dSHong Zhang       PetscBool assembled;
166*33f72c5dSHong Zhang       ierr = MatZeroEntries(Amat);CHKERRQ(ierr);
167*33f72c5dSHong Zhang       ierr = MatAssembled(Amat,&assembled);CHKERRQ(ierr);
168*33f72c5dSHong Zhang       if (!assembled) {
169*33f72c5dSHong Zhang         ierr = MatAssemblyBegin(Amat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
170*33f72c5dSHong Zhang         ierr = MatAssemblyEnd(Amat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
171*33f72c5dSHong Zhang       }
172*33f72c5dSHong Zhang     }
173*33f72c5dSHong Zhang   } else {
174*33f72c5dSHong Zhang     if (ts->rhsjacobianp) {
175*33f72c5dSHong Zhang       ierr = TSComputeRHSJacobianP(ts,t,U,ts->Jacprhs);CHKERRQ(ierr);
176*33f72c5dSHong Zhang     }
177*33f72c5dSHong Zhang     if (ts->Jacprhs == Amat) { /* No IJacobian, so we only have the RHS matrix */
178*33f72c5dSHong Zhang       ierr = MatScale(Amat,-1);CHKERRQ(ierr);
179*33f72c5dSHong Zhang     } else if (ts->Jacprhs) { /* Both IJacobian and RHSJacobian */
180*33f72c5dSHong Zhang       MatStructure axpy = DIFFERENT_NONZERO_PATTERN;
181*33f72c5dSHong Zhang       if (!ts->ijacobianp) { /* No IJacobianp provided, but we have a separate RHS matrix */
182*33f72c5dSHong Zhang         ierr = MatZeroEntries(Amat);CHKERRQ(ierr);
183*33f72c5dSHong Zhang       }
184*33f72c5dSHong Zhang       ierr = MatAXPY(Amat,-1,ts->Jacprhs,axpy);CHKERRQ(ierr);
185*33f72c5dSHong Zhang     }
186*33f72c5dSHong Zhang   }
187*33f72c5dSHong Zhang   ierr = PetscLogEventEnd(TS_JacobianPEval,ts,U,Amat,0);CHKERRQ(ierr);
188*33f72c5dSHong Zhang   PetscFunctionReturn(0);
189*33f72c5dSHong Zhang }
190*33f72c5dSHong Zhang 
191*33f72c5dSHong Zhang /*@C
192814e85d6SHong Zhang     TSSetCostIntegrand - Sets the routine for evaluating the integral term in one or more cost functions
193814e85d6SHong Zhang 
194814e85d6SHong Zhang     Logically Collective on TS
195814e85d6SHong Zhang 
196814e85d6SHong Zhang     Input Parameters:
197814e85d6SHong Zhang +   ts - the TS context obtained from TSCreate()
198814e85d6SHong Zhang .   numcost - number of gradients to be computed, this is the number of cost functions
199814e85d6SHong Zhang .   costintegral - vector that stores the integral values
200814e85d6SHong Zhang .   rf - routine for evaluating the integrand function
201c9ad9de0SHong Zhang .   drduf - function that computes the gradients of the r's with respect to u
202814e85d6SHong Zhang .   drdpf - function that computes the gradients of the r's with respect to p, can be NULL if parametric sensitivity is not desired (mu=NULL)
203814e85d6SHong Zhang .   fwd - flag indicating whether to evaluate cost integral in the forward run or the adjoint run
204814e85d6SHong Zhang -   ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
205814e85d6SHong Zhang 
206814e85d6SHong Zhang     Calling sequence of rf:
207c9ad9de0SHong Zhang $   PetscErrorCode rf(TS ts,PetscReal t,Vec U,Vec F,void *ctx);
208814e85d6SHong Zhang 
209c9ad9de0SHong Zhang     Calling sequence of drduf:
210c9ad9de0SHong Zhang $   PetscErroCode drduf(TS ts,PetscReal t,Vec U,Vec *dRdU,void *ctx);
211814e85d6SHong Zhang 
212814e85d6SHong Zhang     Calling sequence of drdpf:
213c9ad9de0SHong Zhang $   PetscErroCode drdpf(TS ts,PetscReal t,Vec U,Vec *dRdP,void *ctx);
214814e85d6SHong Zhang 
215814e85d6SHong Zhang     Level: intermediate
216814e85d6SHong Zhang 
217814e85d6SHong Zhang     Notes:
218814e85d6SHong Zhang     For optimization there is usually a single cost function (numcost = 1). For sensitivities there may be multiple cost functions
219814e85d6SHong Zhang 
220814e85d6SHong Zhang .keywords: TS, sensitivity analysis, timestep, set, quadrature, function
221814e85d6SHong Zhang 
222814e85d6SHong Zhang .seealso: TSSetRHSJacobianP(), TSGetCostGradients(), TSSetCostGradients()
223814e85d6SHong Zhang @*/
224814e85d6SHong Zhang PetscErrorCode TSSetCostIntegrand(TS ts,PetscInt numcost,Vec costintegral,PetscErrorCode (*rf)(TS,PetscReal,Vec,Vec,void*),
225c9ad9de0SHong Zhang                                                           PetscErrorCode (*drduf)(TS,PetscReal,Vec,Vec*,void*),
226814e85d6SHong Zhang                                                           PetscErrorCode (*drdpf)(TS,PetscReal,Vec,Vec*,void*),
227814e85d6SHong Zhang                                                           PetscBool fwd,void *ctx)
228814e85d6SHong Zhang {
229814e85d6SHong Zhang   PetscErrorCode ierr;
230814e85d6SHong Zhang 
231814e85d6SHong Zhang   PetscFunctionBegin;
232814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
233814e85d6SHong Zhang   if (costintegral) PetscValidHeaderSpecific(costintegral,VEC_CLASSID,3);
234814e85d6SHong Zhang   if (ts->numcost && ts->numcost!=numcost) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"The number of cost functions (2rd parameter of TSSetCostIntegrand()) is inconsistent with the one set by TSSetCostGradients() or TSForwardSetIntegralGradients()");
235814e85d6SHong Zhang   if (!ts->numcost) ts->numcost=numcost;
236814e85d6SHong Zhang 
237814e85d6SHong Zhang   if (costintegral) {
238814e85d6SHong Zhang     ierr = PetscObjectReference((PetscObject)costintegral);CHKERRQ(ierr);
239814e85d6SHong Zhang     ierr = VecDestroy(&ts->vec_costintegral);CHKERRQ(ierr);
240814e85d6SHong Zhang     ts->vec_costintegral = costintegral;
241814e85d6SHong Zhang   } else {
242814e85d6SHong Zhang     if (!ts->vec_costintegral) { /* Create a seq vec if user does not provide one */
243814e85d6SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,numcost,&ts->vec_costintegral);CHKERRQ(ierr);
244814e85d6SHong Zhang     } else {
245814e85d6SHong Zhang       ierr = VecSet(ts->vec_costintegral,0.0);CHKERRQ(ierr);
246814e85d6SHong Zhang     }
247814e85d6SHong Zhang   }
248814e85d6SHong Zhang   if (!ts->vec_costintegrand) {
249814e85d6SHong Zhang     ierr = VecDuplicate(ts->vec_costintegral,&ts->vec_costintegrand);CHKERRQ(ierr);
250814e85d6SHong Zhang   } else {
251814e85d6SHong Zhang     ierr = VecSet(ts->vec_costintegrand,0.0);CHKERRQ(ierr);
252814e85d6SHong Zhang   }
253814e85d6SHong Zhang   ts->costintegralfwd  = fwd; /* Evaluate the cost integral in forward run if fwd is true */
254814e85d6SHong Zhang   ts->costintegrand    = rf;
255814e85d6SHong Zhang   ts->costintegrandctx = ctx;
256c9ad9de0SHong Zhang   ts->drdufunction     = drduf;
257814e85d6SHong Zhang   ts->drdpfunction     = drdpf;
258814e85d6SHong Zhang   PetscFunctionReturn(0);
259814e85d6SHong Zhang }
260814e85d6SHong Zhang 
261b5b4017aSHong Zhang /*@C
262814e85d6SHong Zhang    TSGetCostIntegral - Returns the values of the integral term in the cost functions.
263814e85d6SHong Zhang    It is valid to call the routine after a backward run.
264814e85d6SHong Zhang 
265814e85d6SHong Zhang    Not Collective
266814e85d6SHong Zhang 
267814e85d6SHong Zhang    Input Parameter:
268814e85d6SHong Zhang .  ts - the TS context obtained from TSCreate()
269814e85d6SHong Zhang 
270814e85d6SHong Zhang    Output Parameter:
271814e85d6SHong Zhang .  v - the vector containing the integrals for each cost function
272814e85d6SHong Zhang 
273814e85d6SHong Zhang    Level: intermediate
274814e85d6SHong Zhang 
275814e85d6SHong Zhang .seealso: TSSetCostIntegrand()
276814e85d6SHong Zhang 
277814e85d6SHong Zhang .keywords: TS, sensitivity analysis
278814e85d6SHong Zhang @*/
279814e85d6SHong Zhang PetscErrorCode  TSGetCostIntegral(TS ts,Vec *v)
280814e85d6SHong Zhang {
281814e85d6SHong Zhang   PetscFunctionBegin;
282814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
283814e85d6SHong Zhang   PetscValidPointer(v,2);
284814e85d6SHong Zhang   *v = ts->vec_costintegral;
285814e85d6SHong Zhang   PetscFunctionReturn(0);
286814e85d6SHong Zhang }
287814e85d6SHong Zhang 
288b5b4017aSHong Zhang /*@C
289814e85d6SHong Zhang    TSComputeCostIntegrand - Evaluates the integral function in the cost functions.
290814e85d6SHong Zhang 
291814e85d6SHong Zhang    Input Parameters:
292814e85d6SHong Zhang +  ts - the TS context
293814e85d6SHong Zhang .  t - current time
294c9ad9de0SHong Zhang -  U - state vector, i.e. current solution
295814e85d6SHong Zhang 
296814e85d6SHong Zhang    Output Parameter:
297c9ad9de0SHong Zhang .  Q - vector of size numcost to hold the outputs
298814e85d6SHong Zhang 
299814e85d6SHong Zhang    Note:
300814e85d6SHong Zhang    Most users should not need to explicitly call this routine, as it
301814e85d6SHong Zhang    is used internally within the sensitivity analysis context.
302814e85d6SHong Zhang 
303814e85d6SHong Zhang    Level: developer
304814e85d6SHong Zhang 
305814e85d6SHong Zhang .keywords: TS, compute
306814e85d6SHong Zhang 
307814e85d6SHong Zhang .seealso: TSSetCostIntegrand()
308814e85d6SHong Zhang @*/
309c9ad9de0SHong Zhang PetscErrorCode TSComputeCostIntegrand(TS ts,PetscReal t,Vec U,Vec Q)
310814e85d6SHong Zhang {
311814e85d6SHong Zhang   PetscErrorCode ierr;
312814e85d6SHong Zhang 
313814e85d6SHong Zhang   PetscFunctionBegin;
314814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
315c9ad9de0SHong Zhang   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
316c9ad9de0SHong Zhang   PetscValidHeaderSpecific(Q,VEC_CLASSID,4);
317814e85d6SHong Zhang 
318c9ad9de0SHong Zhang   ierr = PetscLogEventBegin(TS_FunctionEval,ts,U,Q,0);CHKERRQ(ierr);
319814e85d6SHong Zhang   if (ts->costintegrand) {
320814e85d6SHong Zhang     PetscStackPush("TS user integrand in the cost function");
321c9ad9de0SHong Zhang     ierr = (*ts->costintegrand)(ts,t,U,Q,ts->costintegrandctx);CHKERRQ(ierr);
322814e85d6SHong Zhang     PetscStackPop;
323814e85d6SHong Zhang   } else {
324c9ad9de0SHong Zhang     ierr = VecZeroEntries(Q);CHKERRQ(ierr);
325814e85d6SHong Zhang   }
326814e85d6SHong Zhang 
327c9ad9de0SHong Zhang   ierr = PetscLogEventEnd(TS_FunctionEval,ts,U,Q,0);CHKERRQ(ierr);
328814e85d6SHong Zhang   PetscFunctionReturn(0);
329814e85d6SHong Zhang }
330814e85d6SHong Zhang 
331b5b4017aSHong Zhang /*@C
332c9ad9de0SHong Zhang   TSComputeDRDUFunction - Runs the user-defined DRDU function.
333814e85d6SHong Zhang 
334814e85d6SHong Zhang   Collective on TS
335814e85d6SHong Zhang 
336814e85d6SHong Zhang   Input Parameters:
337c9ad9de0SHong Zhang + ts - the TS context obtained from TSCreate()
338c9ad9de0SHong Zhang . t - current time
339c9ad9de0SHong Zhang - U - stata vector
340c9ad9de0SHong Zhang 
341c9ad9de0SHong Zhang   Output Parameters:
342b5b4017aSHong Zhang . DRDU - vector array to hold the outputs
343814e85d6SHong Zhang 
344814e85d6SHong Zhang   Notes:
345c9ad9de0SHong Zhang   TSComputeDRDUFunction() is typically used for sensitivity implementation,
346814e85d6SHong Zhang   so most users would not generally call this routine themselves.
347814e85d6SHong Zhang 
348814e85d6SHong Zhang   Level: developer
349814e85d6SHong Zhang 
350814e85d6SHong Zhang .keywords: TS, sensitivity
351814e85d6SHong Zhang .seealso: TSSetCostIntegrand()
352814e85d6SHong Zhang @*/
353c9ad9de0SHong Zhang PetscErrorCode TSComputeDRDUFunction(TS ts,PetscReal t,Vec U,Vec *DRDU)
354814e85d6SHong Zhang {
355814e85d6SHong Zhang   PetscErrorCode ierr;
356814e85d6SHong Zhang 
357814e85d6SHong Zhang   PetscFunctionBegin;
358*33f72c5dSHong Zhang   if (!DRDU) PetscFunctionReturn(0);
359814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
360c9ad9de0SHong Zhang   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
361814e85d6SHong Zhang 
362c9ad9de0SHong Zhang   PetscStackPush("TS user DRDU function for sensitivity analysis");
363c9ad9de0SHong Zhang   ierr = (*ts->drdufunction)(ts,t,U,DRDU,ts->costintegrandctx); CHKERRQ(ierr);
364814e85d6SHong Zhang   PetscStackPop;
365814e85d6SHong Zhang   PetscFunctionReturn(0);
366814e85d6SHong Zhang }
367814e85d6SHong Zhang 
368b5b4017aSHong Zhang /*@C
369814e85d6SHong Zhang   TSComputeDRDPFunction - Runs the user-defined DRDP function.
370814e85d6SHong Zhang 
371814e85d6SHong Zhang   Collective on TS
372814e85d6SHong Zhang 
373814e85d6SHong Zhang   Input Parameters:
374c9ad9de0SHong Zhang + ts - the TS context obtained from TSCreate()
375c9ad9de0SHong Zhang . t - current time
376c9ad9de0SHong Zhang - U - stata vector
377c9ad9de0SHong Zhang 
378c9ad9de0SHong Zhang   Output Parameters:
379b5b4017aSHong Zhang . DRDP - vector array to hold the outputs
380814e85d6SHong Zhang 
381814e85d6SHong Zhang   Notes:
382c9ad9de0SHong Zhang   TSComputeDRDPFunction() is typically used for sensitivity implementation,
383814e85d6SHong Zhang   so most users would not generally call this routine themselves.
384814e85d6SHong Zhang 
385814e85d6SHong Zhang   Level: developer
386814e85d6SHong Zhang 
387814e85d6SHong Zhang .keywords: TS, sensitivity
388814e85d6SHong Zhang .seealso: TSSetCostIntegrand()
389814e85d6SHong Zhang @*/
390c9ad9de0SHong Zhang PetscErrorCode TSComputeDRDPFunction(TS ts,PetscReal t,Vec U,Vec *DRDP)
391814e85d6SHong Zhang {
392814e85d6SHong Zhang   PetscErrorCode ierr;
393814e85d6SHong Zhang 
394814e85d6SHong Zhang   PetscFunctionBegin;
395*33f72c5dSHong Zhang   if (!DRDP) PetscFunctionReturn(0);
396814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
397c9ad9de0SHong Zhang   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
398814e85d6SHong Zhang 
399814e85d6SHong Zhang   PetscStackPush("TS user DRDP function for sensitivity analysis");
400c9ad9de0SHong Zhang   ierr = (*ts->drdpfunction)(ts,t,U,DRDP,ts->costintegrandctx); CHKERRQ(ierr);
401814e85d6SHong Zhang   PetscStackPop;
402814e85d6SHong Zhang   PetscFunctionReturn(0);
403814e85d6SHong Zhang }
404814e85d6SHong Zhang 
405b5b4017aSHong Zhang /*@C
406b5b4017aSHong 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.
407b5b4017aSHong Zhang 
408b5b4017aSHong Zhang   Logically Collective on TS
409b5b4017aSHong Zhang 
410b5b4017aSHong Zhang   Input Parameters:
411b5b4017aSHong Zhang + ts - TS context obtained from TSCreate()
412b5b4017aSHong Zhang . ihp1 - an array of vectors storing the result of vector-Hessian-vector product for F_UU
413b5b4017aSHong Zhang . hessianproductfunc1 - vector-Hessian-vector product function for F_UU
414b5b4017aSHong Zhang . ihp2 - an array of vectors storing the result of vector-Hessian-vector product for F_UP
415b5b4017aSHong Zhang . hessianproductfunc2 - vector-Hessian-vector product function for F_UP
416b5b4017aSHong Zhang . ihp3 - an array of vectors storing the result of vector-Hessian-vector product for F_PU
417b5b4017aSHong Zhang . hessianproductfunc3 - vector-Hessian-vector product function for F_PU
418b5b4017aSHong Zhang . ihp4 - an array of vectors storing the result of vector-Hessian-vector product for F_PP
419b5b4017aSHong Zhang . hessianproductfunc4 - vector-Hessian-vector product function for F_PP
420b5b4017aSHong Zhang 
421b5b4017aSHong Zhang   Calling sequence of ihessianproductfunc:
422b5b4017aSHong Zhang $ ihessianproductfunc (TS ts,PetscReal t,Vec U,Vec Vl,Vec Vr,Vec VHV,void *ctx);
423b5b4017aSHong Zhang +   t - current timestep
424b5b4017aSHong Zhang .   U - input vector (current ODE solution)
425b5b4017aSHong Zhang .   Vl - input vector to be left-multiplied with the Hessian
426b5b4017aSHong Zhang .   Vr - input vector to be right-multiplied with the Hessian
427b5b4017aSHong Zhang .   VHV - output vector for vector-Hessian-vector product
428b5b4017aSHong Zhang -   ctx - [optional] user-defined function context
429b5b4017aSHong Zhang 
430b5b4017aSHong Zhang   Level: intermediate
431b5b4017aSHong Zhang 
432b5b4017aSHong Zhang Note: The first Hessian function and the working array are required.
433b5b4017aSHong Zhang 
434b5b4017aSHong Zhang .keywords: TS, sensitivity
435b5b4017aSHong Zhang 
436b5b4017aSHong Zhang .seealso:
437b5b4017aSHong Zhang @*/
438b5b4017aSHong Zhang PetscErrorCode TSSetIHessianProduct(TS ts,Vec *ihp1,PetscErrorCode (*ihessianproductfunc1)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*),
439b5b4017aSHong Zhang                                           Vec *ihp2,PetscErrorCode (*ihessianproductfunc2)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*),
440b5b4017aSHong Zhang                                           Vec *ihp3,PetscErrorCode (*ihessianproductfunc3)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*),
441b5b4017aSHong Zhang                                           Vec *ihp4,PetscErrorCode (*ihessianproductfunc4)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*),
442b5b4017aSHong Zhang                                     void *ctx)
443b5b4017aSHong Zhang {
444b5b4017aSHong Zhang   PetscFunctionBegin;
445b5b4017aSHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
446b5b4017aSHong Zhang   PetscValidPointer(ihp1,2);
447b5b4017aSHong Zhang 
448b5b4017aSHong Zhang   ts->ihessianproductctx = ctx;
449b5b4017aSHong Zhang   if (ihp1) ts->vecs_fuu = ihp1;
450b5b4017aSHong Zhang   if (ihp2) ts->vecs_fup = ihp2;
451b5b4017aSHong Zhang   if (ihp3) ts->vecs_fpu = ihp3;
452b5b4017aSHong Zhang   if (ihp4) ts->vecs_fpp = ihp4;
453b5b4017aSHong Zhang   ts->ihessianproduct_fuu = ihessianproductfunc1;
454b5b4017aSHong Zhang   ts->ihessianproduct_fup = ihessianproductfunc2;
455b5b4017aSHong Zhang   ts->ihessianproduct_fpu = ihessianproductfunc3;
456b5b4017aSHong Zhang   ts->ihessianproduct_fpp = ihessianproductfunc4;
457b5b4017aSHong Zhang   PetscFunctionReturn(0);
458b5b4017aSHong Zhang }
459b5b4017aSHong Zhang 
460b5b4017aSHong Zhang /*@C
461b5b4017aSHong Zhang   TSComputeIHessianProductFunction1 - Runs the user-defined vector-Hessian-vector product function for Fuu.
462b5b4017aSHong Zhang 
463b5b4017aSHong Zhang   Collective on TS
464b5b4017aSHong Zhang 
465b5b4017aSHong Zhang   Input Parameters:
466b5b4017aSHong Zhang . ts   - The TS context obtained from TSCreate()
467b5b4017aSHong Zhang 
468b5b4017aSHong Zhang   Notes:
469b5b4017aSHong Zhang   TSComputeIHessianProductFunction1() is typically used for sensitivity implementation,
470b5b4017aSHong Zhang   so most users would not generally call this routine themselves.
471b5b4017aSHong Zhang 
472b5b4017aSHong Zhang   Level: developer
473b5b4017aSHong Zhang 
474b5b4017aSHong Zhang .keywords: TS, sensitivity
475b5b4017aSHong Zhang 
476b5b4017aSHong Zhang .seealso: TSSetIHessianProduct()
477b5b4017aSHong Zhang @*/
478b5b4017aSHong Zhang PetscErrorCode TSComputeIHessianProductFunction1(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV)
479b5b4017aSHong Zhang {
480b5b4017aSHong Zhang   PetscErrorCode ierr;
481b5b4017aSHong Zhang 
482b5b4017aSHong Zhang   PetscFunctionBegin;
483*33f72c5dSHong Zhang   if (!VHV) PetscFunctionReturn(0);
484b5b4017aSHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
485b5b4017aSHong Zhang   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
486b5b4017aSHong Zhang 
487b5b4017aSHong Zhang   PetscStackPush("TS user IHessianProduct function 1 for sensitivity analysis");
488b5b4017aSHong Zhang   ierr = (*ts->ihessianproduct_fuu)(ts,t,U,Vl,Vr,VHV,ts->ihessianproductctx);CHKERRQ(ierr);
489b5b4017aSHong Zhang   PetscStackPop;
490b5b4017aSHong Zhang   PetscFunctionReturn(0);
491b5b4017aSHong Zhang }
492b5b4017aSHong Zhang 
493b5b4017aSHong Zhang /*@C
494b5b4017aSHong Zhang   TSComputeIHessianProductFunction2 - Runs the user-defined vector-Hessian-vector product function for Fup.
495b5b4017aSHong Zhang 
496b5b4017aSHong Zhang   Collective on TS
497b5b4017aSHong Zhang 
498b5b4017aSHong Zhang   Input Parameters:
499b5b4017aSHong Zhang . ts   - The TS context obtained from TSCreate()
500b5b4017aSHong Zhang 
501b5b4017aSHong Zhang   Notes:
502b5b4017aSHong Zhang   TSComputeIHessianProductFunction2() is typically used for sensitivity implementation,
503b5b4017aSHong Zhang   so most users would not generally call this routine themselves.
504b5b4017aSHong Zhang 
505b5b4017aSHong Zhang   Level: developer
506b5b4017aSHong Zhang 
507b5b4017aSHong Zhang .keywords: TS, sensitivity
508b5b4017aSHong Zhang 
509b5b4017aSHong Zhang .seealso: TSSetIHessianProduct()
510b5b4017aSHong Zhang @*/
511b5b4017aSHong Zhang PetscErrorCode TSComputeIHessianProductFunction2(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV)
512b5b4017aSHong Zhang {
513b5b4017aSHong Zhang   PetscErrorCode ierr;
514b5b4017aSHong Zhang 
515b5b4017aSHong Zhang   PetscFunctionBegin;
516*33f72c5dSHong Zhang   if (!VHV) PetscFunctionReturn(0);
517b5b4017aSHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
518b5b4017aSHong Zhang   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
519b5b4017aSHong Zhang 
520b5b4017aSHong Zhang   PetscStackPush("TS user IHessianProduct function 2 for sensitivity analysis");
521b5b4017aSHong Zhang   ierr = (*ts->ihessianproduct_fup)(ts,t,U,Vl,Vr,VHV,ts->ihessianproductctx);CHKERRQ(ierr);
522b5b4017aSHong Zhang   PetscStackPop;
523b5b4017aSHong Zhang   PetscFunctionReturn(0);
524b5b4017aSHong Zhang }
525b5b4017aSHong Zhang 
526b5b4017aSHong Zhang /*@C
527b5b4017aSHong Zhang   TSComputeIHessianProductFunction3 - Runs the user-defined vector-Hessian-vector product function for Fpu.
528b5b4017aSHong Zhang 
529b5b4017aSHong Zhang   Collective on TS
530b5b4017aSHong Zhang 
531b5b4017aSHong Zhang   Input Parameters:
532b5b4017aSHong Zhang . ts   - The TS context obtained from TSCreate()
533b5b4017aSHong Zhang 
534b5b4017aSHong Zhang   Notes:
535b5b4017aSHong Zhang   TSComputeIHessianProductFunction3() is typically used for sensitivity implementation,
536b5b4017aSHong Zhang   so most users would not generally call this routine themselves.
537b5b4017aSHong Zhang 
538b5b4017aSHong Zhang   Level: developer
539b5b4017aSHong Zhang 
540b5b4017aSHong Zhang .keywords: TS, sensitivity
541b5b4017aSHong Zhang 
542b5b4017aSHong Zhang .seealso: TSSetIHessianProduct()
543b5b4017aSHong Zhang @*/
544b5b4017aSHong Zhang PetscErrorCode TSComputeIHessianProductFunction3(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV)
545b5b4017aSHong Zhang {
546b5b4017aSHong Zhang   PetscErrorCode ierr;
547b5b4017aSHong Zhang 
548b5b4017aSHong Zhang   PetscFunctionBegin;
549*33f72c5dSHong Zhang   if (!VHV) PetscFunctionReturn(0);
550b5b4017aSHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
551b5b4017aSHong Zhang   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
552b5b4017aSHong Zhang 
553b5b4017aSHong Zhang   PetscStackPush("TS user IHessianProduct function 3 for sensitivity analysis");
554b5b4017aSHong Zhang   ierr = (*ts->ihessianproduct_fpu)(ts,t,U,Vl,Vr,VHV,ts->ihessianproductctx);CHKERRQ(ierr);
555b5b4017aSHong Zhang   PetscStackPop;
556b5b4017aSHong Zhang   PetscFunctionReturn(0);
557b5b4017aSHong Zhang }
558b5b4017aSHong Zhang 
559b5b4017aSHong Zhang /*@C
560b5b4017aSHong Zhang   TSComputeIHessianProductFunction4 - Runs the user-defined vector-Hessian-vector product function for Fpp.
561b5b4017aSHong Zhang 
562b5b4017aSHong Zhang   Collective on TS
563b5b4017aSHong Zhang 
564b5b4017aSHong Zhang   Input Parameters:
565b5b4017aSHong Zhang . ts   - The TS context obtained from TSCreate()
566b5b4017aSHong Zhang 
567b5b4017aSHong Zhang   Notes:
568b5b4017aSHong Zhang   TSComputeIHessianProductFunction4() is typically used for sensitivity implementation,
569b5b4017aSHong Zhang   so most users would not generally call this routine themselves.
570b5b4017aSHong Zhang 
571b5b4017aSHong Zhang   Level: developer
572b5b4017aSHong Zhang 
573b5b4017aSHong Zhang .keywords: TS, sensitivity
574b5b4017aSHong Zhang 
575b5b4017aSHong Zhang .seealso: TSSetIHessianProduct()
576b5b4017aSHong Zhang @*/
577b5b4017aSHong Zhang PetscErrorCode TSComputeIHessianProductFunction4(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV)
578b5b4017aSHong Zhang {
579b5b4017aSHong Zhang   PetscErrorCode ierr;
580b5b4017aSHong Zhang 
581b5b4017aSHong Zhang   PetscFunctionBegin;
582*33f72c5dSHong Zhang   if (!VHV) PetscFunctionReturn(0);
583b5b4017aSHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
584b5b4017aSHong Zhang   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
585b5b4017aSHong Zhang 
586b5b4017aSHong Zhang   PetscStackPush("TS user IHessianProduct function 3 for sensitivity analysis");
587b5b4017aSHong Zhang   ierr = (*ts->ihessianproduct_fpp)(ts,t,U,Vl,Vr,VHV,ts->ihessianproductctx);CHKERRQ(ierr);
588b5b4017aSHong Zhang   PetscStackPop;
589b5b4017aSHong Zhang   PetscFunctionReturn(0);
590b5b4017aSHong Zhang }
591b5b4017aSHong Zhang 
592814e85d6SHong Zhang /* --------------------------- Adjoint sensitivity ---------------------------*/
593814e85d6SHong Zhang 
594814e85d6SHong Zhang /*@
595814e85d6SHong 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
596814e85d6SHong Zhang       for use by the TSAdjoint routines.
597814e85d6SHong Zhang 
598814e85d6SHong Zhang    Logically Collective on TS and Vec
599814e85d6SHong Zhang 
600814e85d6SHong Zhang    Input Parameters:
601814e85d6SHong Zhang +  ts - the TS context obtained from TSCreate()
602814e85d6SHong 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
603814e85d6SHong Zhang -  mu - gradients with respect to the parameters, the number of entries in these vectors is the same as the number of parameters
604814e85d6SHong Zhang 
605814e85d6SHong Zhang    Level: beginner
606814e85d6SHong Zhang 
607814e85d6SHong Zhang    Notes:
608814e85d6SHong Zhang     the entries in these vectors must be correctly initialized with the values lamda_i = df/dy|finaltime  mu_i = df/dp|finaltime
609814e85d6SHong Zhang 
610814e85d6SHong Zhang    After TSAdjointSolve() is called the lamba and the mu contain the computed sensitivities
611814e85d6SHong Zhang 
612814e85d6SHong Zhang .keywords: TS, timestep, set, sensitivity, initial values
613b5b4017aSHong Zhang 
614b5b4017aSHong Zhang .seealso TSGetCostGradients()
615814e85d6SHong Zhang @*/
616814e85d6SHong Zhang PetscErrorCode TSSetCostGradients(TS ts,PetscInt numcost,Vec *lambda,Vec *mu)
617814e85d6SHong Zhang {
618814e85d6SHong Zhang   PetscFunctionBegin;
619814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
620814e85d6SHong Zhang   PetscValidPointer(lambda,2);
621814e85d6SHong Zhang   ts->vecs_sensi  = lambda;
622814e85d6SHong Zhang   ts->vecs_sensip = mu;
623814e85d6SHong Zhang   if (ts->numcost && ts->numcost!=numcost) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"The number of cost functions (2rd parameter of TSSetCostIntegrand()) is inconsistent with the one set by TSSetCostIntegrand");
624814e85d6SHong Zhang   ts->numcost  = numcost;
625814e85d6SHong Zhang   PetscFunctionReturn(0);
626814e85d6SHong Zhang }
627814e85d6SHong Zhang 
628814e85d6SHong Zhang /*@
629814e85d6SHong Zhang    TSGetCostGradients - Returns the gradients from the TSAdjointSolve()
630814e85d6SHong Zhang 
631814e85d6SHong Zhang    Not Collective, but Vec returned is parallel if TS is parallel
632814e85d6SHong Zhang 
633814e85d6SHong Zhang    Input Parameter:
634814e85d6SHong Zhang .  ts - the TS context obtained from TSCreate()
635814e85d6SHong Zhang 
636814e85d6SHong Zhang    Output Parameter:
637814e85d6SHong Zhang +  lambda - vectors containing the gradients of the cost functions with respect to the ODE/DAE solution variables
638814e85d6SHong Zhang -  mu - vectors containing the gradients of the cost functions with respect to the problem parameters
639814e85d6SHong Zhang 
640814e85d6SHong Zhang    Level: intermediate
641814e85d6SHong Zhang 
642814e85d6SHong Zhang .keywords: TS, timestep, get, sensitivity
643b5b4017aSHong Zhang 
644b5b4017aSHong Zhang .seealso: TSSetCostGradients()
645814e85d6SHong Zhang @*/
646814e85d6SHong Zhang PetscErrorCode TSGetCostGradients(TS ts,PetscInt *numcost,Vec **lambda,Vec **mu)
647814e85d6SHong Zhang {
648814e85d6SHong Zhang   PetscFunctionBegin;
649814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
650814e85d6SHong Zhang   if (numcost) *numcost = ts->numcost;
651814e85d6SHong Zhang   if (lambda)  *lambda  = ts->vecs_sensi;
652814e85d6SHong Zhang   if (mu)      *mu      = ts->vecs_sensip;
653814e85d6SHong Zhang   PetscFunctionReturn(0);
654814e85d6SHong Zhang }
655814e85d6SHong Zhang 
656814e85d6SHong Zhang /*@
657b5b4017aSHong 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
658b5b4017aSHong Zhang       for use by the TSAdjoint routines.
659b5b4017aSHong Zhang 
660b5b4017aSHong Zhang    Logically Collective on TS and Vec
661b5b4017aSHong Zhang 
662b5b4017aSHong Zhang    Input Parameters:
663b5b4017aSHong Zhang +  ts - the TS context obtained from TSCreate()
664b5b4017aSHong Zhang .  numcost - number of cost functions
665b5b4017aSHong 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
666b5b4017aSHong 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
667b5b4017aSHong Zhang -  dir - the direction vector that are multiplied with the Hessian of the cost functions
668b5b4017aSHong Zhang 
669b5b4017aSHong Zhang    Level: beginner
670b5b4017aSHong Zhang 
671b5b4017aSHong Zhang    Notes: Hessian of the cost function is completely different from Hessian of the ODE/DAE system
672b5b4017aSHong Zhang 
673b5b4017aSHong Zhang    For second-order adjoint, one needs to call this function and then TSAdjointInitializeForward() before TSSolve().
674b5b4017aSHong Zhang 
675b5b4017aSHong Zhang    After TSAdjointSolve() is called, the lamba2 and the mu2 will contain the computed second-order adjoint sensitivities, and can be used to produce Hessian-vector product (not the full Hessian matrix). Users must provide a direction vector; it is usually generated by an optimization solver.
676b5b4017aSHong Zhang 
6773fca9621SHong Zhang    Passing NULL for lambda2 disables the second-order calculation.
678b5b4017aSHong Zhang .keywords: TS, sensitivity, second-order adjoint
679b5b4017aSHong Zhang 
680b5b4017aSHong Zhang .seealso: TSAdjointInitializeForward()
681b5b4017aSHong Zhang @*/
682b5b4017aSHong Zhang PetscErrorCode TSSetCostHessianProducts(TS ts,PetscInt numcost,Vec *lambda2,Vec *mu2,Vec dir)
683b5b4017aSHong Zhang {
684b5b4017aSHong Zhang   PetscFunctionBegin;
685b5b4017aSHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
686b5b4017aSHong Zhang   if (ts->numcost && ts->numcost!=numcost) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"The number of cost functions (2rd parameter of TSSetCostIntegrand()) is inconsistent with the one set by TSSetCostIntegrand");
687b5b4017aSHong Zhang   ts->numcost       = numcost;
688b5b4017aSHong Zhang   ts->vecs_sensi2   = lambda2;
689*33f72c5dSHong Zhang   ts->vecs_sensi2p  = mu2;
690b5b4017aSHong Zhang   ts->vec_dir       = dir;
691b5b4017aSHong Zhang   PetscFunctionReturn(0);
692b5b4017aSHong Zhang }
693b5b4017aSHong Zhang 
694b5b4017aSHong Zhang /*@
695b5b4017aSHong Zhang    TSGetCostHessianProducts - Returns the gradients from the TSAdjointSolve()
696b5b4017aSHong Zhang 
697b5b4017aSHong Zhang    Not Collective, but Vec returned is parallel if TS is parallel
698b5b4017aSHong Zhang 
699b5b4017aSHong Zhang    Input Parameter:
700b5b4017aSHong Zhang .  ts - the TS context obtained from TSCreate()
701b5b4017aSHong Zhang 
702b5b4017aSHong Zhang    Output Parameter:
703b5b4017aSHong Zhang +  numcost - number of cost functions
704b5b4017aSHong 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
705b5b4017aSHong 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
706b5b4017aSHong Zhang -  dir - the direction vector that are multiplied with the Hessian of the cost functions
707b5b4017aSHong Zhang 
708b5b4017aSHong Zhang    Level: intermediate
709b5b4017aSHong Zhang 
710b5b4017aSHong Zhang .keywords: TS, sensitivity, second-order adjoint
711b5b4017aSHong Zhang 
712b5b4017aSHong Zhang .seealso: TSSetCostHessianProducts()
713b5b4017aSHong Zhang @*/
714b5b4017aSHong Zhang PetscErrorCode TSGetCostHessianProducts(TS ts,PetscInt *numcost,Vec **lambda2,Vec **mu2, Vec *dir)
715b5b4017aSHong Zhang {
716b5b4017aSHong Zhang   PetscFunctionBegin;
717b5b4017aSHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
718b5b4017aSHong Zhang   if (numcost) *numcost = ts->numcost;
719b5b4017aSHong Zhang   if (lambda2) *lambda2 = ts->vecs_sensi2;
720*33f72c5dSHong Zhang   if (mu2)     *mu2     = ts->vecs_sensi2p;
721b5b4017aSHong Zhang   if (dir)     *dir     = ts->vec_dir;
722b5b4017aSHong Zhang   PetscFunctionReturn(0);
723b5b4017aSHong Zhang }
724b5b4017aSHong Zhang 
725b5b4017aSHong Zhang /*@
726b5b4017aSHong Zhang   TSAdjointInitializeForward - Trigger the tangent linear solver and initialize the forward sensitivities
727b5b4017aSHong Zhang 
728b5b4017aSHong Zhang   Logically Collective on TS and Mat
729b5b4017aSHong Zhang 
730b5b4017aSHong Zhang   Input Parameters:
731b5b4017aSHong Zhang +  ts - the TS context obtained from TSCreate()
732b5b4017aSHong Zhang -  didp - the derivative of initial values w.r.t. parameters
733b5b4017aSHong Zhang 
734b5b4017aSHong Zhang   Level: intermediate
735b5b4017aSHong Zhang 
736b5b4017aSHong Zhang   Notes: When computing sensitivies w.r.t. initial condition, set didp to NULL so that the solver will take it as an identity matrix mathematically.
737b5b4017aSHong Zhang 
738b5b4017aSHong Zhang .keywords: TS, sensitivity, second-order adjoint
739b5b4017aSHong Zhang 
740b5b4017aSHong Zhang .seealso: TSSetCostHessianProducts()
741b5b4017aSHong Zhang @*/
742b5b4017aSHong Zhang PetscErrorCode TSAdjointInitializeForward(TS ts,Mat didp)
743b5b4017aSHong Zhang {
744*33f72c5dSHong Zhang   Mat            A;
745*33f72c5dSHong Zhang   Vec            sp;
746*33f72c5dSHong Zhang   PetscScalar    *xarr;
747*33f72c5dSHong Zhang   PetscInt       lsize;
748b5b4017aSHong Zhang   PetscErrorCode ierr;
749b5b4017aSHong Zhang 
750b5b4017aSHong Zhang   PetscFunctionBegin;
751b5b4017aSHong Zhang   ts->forward_solve = PETSC_TRUE; /* turn on tangent linear mode */
752b5b4017aSHong Zhang   if (!ts->vecs_sensi2) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must call TSSetCostHessianProducts() first");
753*33f72c5dSHong Zhang   if (!ts->vec_dir) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Directional vector is missing. Call TSSetCostHessianProducts() to set it.");
754*33f72c5dSHong Zhang   /* create a single-column dense matrix */
755*33f72c5dSHong Zhang   ierr = VecGetLocalSize(ts->vec_sol,&lsize);CHKERRQ(ierr);
756*33f72c5dSHong Zhang   ierr = MatCreateDense(PETSC_COMM_WORLD,lsize,PETSC_DECIDE,PETSC_DECIDE,1,NULL,&A);CHKERRQ(ierr);
757*33f72c5dSHong Zhang 
758*33f72c5dSHong Zhang   ierr = VecDuplicate(ts->vec_sol,&sp);CHKERRQ(ierr);
759*33f72c5dSHong Zhang   ierr = MatDenseGetColumn(A,0,&xarr);CHKERRQ(ierr);
760*33f72c5dSHong Zhang   ierr = VecPlaceArray(sp,xarr);CHKERRQ(ierr);
761*33f72c5dSHong Zhang   if (ts->vecs_sensi2p) { /* TLM variable initialized as 2*dIdP*dir */
762*33f72c5dSHong Zhang     if (didp) {
763*33f72c5dSHong Zhang       ierr = MatMult(didp,ts->vec_dir,sp);CHKERRQ(ierr);
764*33f72c5dSHong Zhang       ierr = VecScale(sp,2.);CHKERRQ(ierr);
765*33f72c5dSHong Zhang     } else {
766*33f72c5dSHong Zhang       ierr = VecZeroEntries(sp);CHKERRQ(ierr);
767*33f72c5dSHong Zhang     }
768*33f72c5dSHong Zhang   } else { /* TLM variable initialized as dir */
769*33f72c5dSHong Zhang     ierr = VecCopy(ts->vec_dir,sp);CHKERRQ(ierr);
770*33f72c5dSHong Zhang   }
771*33f72c5dSHong Zhang   ierr = VecResetArray(sp);CHKERRQ(ierr);
772*33f72c5dSHong Zhang   ierr = MatDenseRestoreColumn(A,&xarr);CHKERRQ(ierr);
773*33f72c5dSHong Zhang   ierr = VecDestroy(&sp);CHKERRQ(ierr);
774*33f72c5dSHong Zhang 
775*33f72c5dSHong Zhang   ierr = TSForwardSetInitialSensitivities(ts,A);CHKERRQ(ierr); /* if didp is NULL, identity matrix is assumed */
776*33f72c5dSHong Zhang 
777*33f72c5dSHong Zhang   ierr = MatDestroy(&A);CHKERRQ(ierr);
778b5b4017aSHong Zhang   PetscFunctionReturn(0);
779b5b4017aSHong Zhang }
780b5b4017aSHong Zhang 
781b5b4017aSHong Zhang /*@
782814e85d6SHong Zhang    TSAdjointSetUp - Sets up the internal data structures for the later use
783814e85d6SHong Zhang    of an adjoint solver
784814e85d6SHong Zhang 
785814e85d6SHong Zhang    Collective on TS
786814e85d6SHong Zhang 
787814e85d6SHong Zhang    Input Parameter:
788814e85d6SHong Zhang .  ts - the TS context obtained from TSCreate()
789814e85d6SHong Zhang 
790814e85d6SHong Zhang    Level: advanced
791814e85d6SHong Zhang 
792814e85d6SHong Zhang .keywords: TS, timestep, setup
793814e85d6SHong Zhang 
794814e85d6SHong Zhang .seealso: TSCreate(), TSAdjointStep(), TSSetCostGradients()
795814e85d6SHong Zhang @*/
796814e85d6SHong Zhang PetscErrorCode TSAdjointSetUp(TS ts)
797814e85d6SHong Zhang {
798814e85d6SHong Zhang   PetscErrorCode ierr;
799814e85d6SHong Zhang 
800814e85d6SHong Zhang   PetscFunctionBegin;
801814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
802814e85d6SHong Zhang   if (ts->adjointsetupcalled) PetscFunctionReturn(0);
803814e85d6SHong Zhang   if (!ts->vecs_sensi) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetCostGradients() first");
804*33f72c5dSHong Zhang   if (ts->vecs_sensip && !ts->Jacp && !ts->Jacprhs) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetRHSJacobianP() or TSSetIJacobianP() first");
805*33f72c5dSHong Zhang 
806*33f72c5dSHong Zhang   if (!ts->Jacp && ts->Jacprhs) ts->Jacp = ts->Jacprhs;
807814e85d6SHong Zhang 
808814e85d6SHong Zhang   if (ts->vec_costintegral) { /* if there is integral in the cost function */
809c9ad9de0SHong Zhang     ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&ts->vecs_drdu);CHKERRQ(ierr);
810814e85d6SHong Zhang     if (ts->vecs_sensip){
811814e85d6SHong Zhang       ierr = VecDuplicateVecs(ts->vecs_sensip[0],ts->numcost,&ts->vecs_drdp);CHKERRQ(ierr);
812814e85d6SHong Zhang     }
813814e85d6SHong Zhang   }
814814e85d6SHong Zhang 
815814e85d6SHong Zhang   if (ts->ops->adjointsetup) {
816814e85d6SHong Zhang     ierr = (*ts->ops->adjointsetup)(ts);CHKERRQ(ierr);
817814e85d6SHong Zhang   }
818814e85d6SHong Zhang   ts->adjointsetupcalled = PETSC_TRUE;
819814e85d6SHong Zhang   PetscFunctionReturn(0);
820814e85d6SHong Zhang }
821814e85d6SHong Zhang 
822814e85d6SHong Zhang /*@
823ece46509SHong Zhang    TSAdjointReset - Resets a TSAdjoint context and removes any allocated Vecs and Mats.
824ece46509SHong Zhang 
825ece46509SHong Zhang    Collective on TS
826ece46509SHong Zhang 
827ece46509SHong Zhang    Input Parameter:
828ece46509SHong Zhang .  ts - the TS context obtained from TSCreate()
829ece46509SHong Zhang 
830ece46509SHong Zhang    Level: beginner
831ece46509SHong Zhang 
832ece46509SHong Zhang .keywords: TS, timestep, reset
833ece46509SHong Zhang 
834ece46509SHong Zhang .seealso: TSCreate(), TSAdjointSetup(), TSADestroy()
835ece46509SHong Zhang @*/
836ece46509SHong Zhang PetscErrorCode TSAdjointReset(TS ts)
837ece46509SHong Zhang {
838ece46509SHong Zhang   PetscErrorCode ierr;
839ece46509SHong Zhang 
840ece46509SHong Zhang   PetscFunctionBegin;
841ece46509SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
842ece46509SHong Zhang   if (ts->ops->adjointreset) {
843ece46509SHong Zhang     ierr = (*ts->ops->adjointreset)(ts);CHKERRQ(ierr);
844ece46509SHong Zhang   }
845cda2db4bSHong Zhang   if (ts->vec_dir) { /* second-order adjoint */
846cda2db4bSHong Zhang     ierr = TSForwardReset(ts);CHKERRQ(ierr);
847cda2db4bSHong Zhang   }
848ece46509SHong Zhang   ts->vecs_sensi         = NULL;
849ece46509SHong Zhang   ts->vecs_sensip        = NULL;
850ece46509SHong Zhang   ts->vecs_sensi2        = NULL;
851*33f72c5dSHong Zhang   ts->vecs_sensi2p       = NULL;
852ece46509SHong Zhang   ts->vec_dir            = NULL;
853ece46509SHong Zhang   ts->adjointsetupcalled = PETSC_FALSE;
854ece46509SHong Zhang   PetscFunctionReturn(0);
855ece46509SHong Zhang }
856ece46509SHong Zhang 
857ece46509SHong Zhang /*@
858814e85d6SHong Zhang    TSAdjointSetSteps - Sets the number of steps the adjoint solver should take backward in time
859814e85d6SHong Zhang 
860814e85d6SHong Zhang    Logically Collective on TS
861814e85d6SHong Zhang 
862814e85d6SHong Zhang    Input Parameters:
863814e85d6SHong Zhang +  ts - the TS context obtained from TSCreate()
864814e85d6SHong Zhang .  steps - number of steps to use
865814e85d6SHong Zhang 
866814e85d6SHong Zhang    Level: intermediate
867814e85d6SHong Zhang 
868814e85d6SHong Zhang    Notes:
869814e85d6SHong Zhang     Normally one does not call this and TSAdjointSolve() integrates back to the original timestep. One can call this
870814e85d6SHong Zhang           so as to integrate back to less than the original timestep
871814e85d6SHong Zhang 
872814e85d6SHong Zhang .keywords: TS, timestep, set, maximum, iterations
873814e85d6SHong Zhang 
874814e85d6SHong Zhang .seealso: TSSetExactFinalTime()
875814e85d6SHong Zhang @*/
876814e85d6SHong Zhang PetscErrorCode TSAdjointSetSteps(TS ts,PetscInt steps)
877814e85d6SHong Zhang {
878814e85d6SHong Zhang   PetscFunctionBegin;
879814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
880814e85d6SHong Zhang   PetscValidLogicalCollectiveInt(ts,steps,2);
881814e85d6SHong Zhang   if (steps < 0) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Cannot step back a negative number of steps");
882814e85d6SHong Zhang   if (steps > ts->steps) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Cannot step back more than the total number of forward steps");
883814e85d6SHong Zhang   ts->adjoint_max_steps = steps;
884814e85d6SHong Zhang   PetscFunctionReturn(0);
885814e85d6SHong Zhang }
886814e85d6SHong Zhang 
887814e85d6SHong Zhang /*@C
888814e85d6SHong Zhang   TSAdjointSetRHSJacobian - Deprecated, use TSSetRHSJacobianP()
889814e85d6SHong Zhang 
890814e85d6SHong Zhang   Level: deprecated
891814e85d6SHong Zhang 
892814e85d6SHong Zhang @*/
893814e85d6SHong Zhang PetscErrorCode TSAdjointSetRHSJacobian(TS ts,Mat Amat,PetscErrorCode (*func)(TS,PetscReal,Vec,Mat,void*),void *ctx)
894814e85d6SHong Zhang {
895814e85d6SHong Zhang   PetscErrorCode ierr;
896814e85d6SHong Zhang 
897814e85d6SHong Zhang   PetscFunctionBegin;
898814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
899814e85d6SHong Zhang   PetscValidHeaderSpecific(Amat,MAT_CLASSID,2);
900814e85d6SHong Zhang 
901814e85d6SHong Zhang   ts->rhsjacobianp    = func;
902814e85d6SHong Zhang   ts->rhsjacobianpctx = ctx;
903814e85d6SHong Zhang   if(Amat) {
904814e85d6SHong Zhang     ierr = PetscObjectReference((PetscObject)Amat);CHKERRQ(ierr);
905814e85d6SHong Zhang     ierr = MatDestroy(&ts->Jacp);CHKERRQ(ierr);
906814e85d6SHong Zhang     ts->Jacp = Amat;
907814e85d6SHong Zhang   }
908814e85d6SHong Zhang   PetscFunctionReturn(0);
909814e85d6SHong Zhang }
910814e85d6SHong Zhang 
911814e85d6SHong Zhang /*@C
912814e85d6SHong Zhang   TSAdjointComputeRHSJacobian - Deprecated, use TSComputeRHSJacobianP()
913814e85d6SHong Zhang 
914814e85d6SHong Zhang   Level: deprecated
915814e85d6SHong Zhang 
916814e85d6SHong Zhang @*/
917c9ad9de0SHong Zhang PetscErrorCode TSAdjointComputeRHSJacobian(TS ts,PetscReal t,Vec U,Mat Amat)
918814e85d6SHong Zhang {
919814e85d6SHong Zhang   PetscErrorCode ierr;
920814e85d6SHong Zhang 
921814e85d6SHong Zhang   PetscFunctionBegin;
922814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
923c9ad9de0SHong Zhang   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
924814e85d6SHong Zhang   PetscValidPointer(Amat,4);
925814e85d6SHong Zhang 
926814e85d6SHong Zhang   PetscStackPush("TS user JacobianP function for sensitivity analysis");
927c9ad9de0SHong Zhang   ierr = (*ts->rhsjacobianp)(ts,t,U,Amat,ts->rhsjacobianpctx); CHKERRQ(ierr);
928814e85d6SHong Zhang   PetscStackPop;
929814e85d6SHong Zhang   PetscFunctionReturn(0);
930814e85d6SHong Zhang }
931814e85d6SHong Zhang 
932814e85d6SHong Zhang /*@
933c9ad9de0SHong Zhang   TSAdjointComputeDRDYFunction - Deprecated, use TSComputeDRDUFunction()
934814e85d6SHong Zhang 
935814e85d6SHong Zhang   Level: deprecated
936814e85d6SHong Zhang 
937814e85d6SHong Zhang @*/
938c9ad9de0SHong Zhang PetscErrorCode TSAdjointComputeDRDYFunction(TS ts,PetscReal t,Vec U,Vec *DRDU)
939814e85d6SHong Zhang {
940814e85d6SHong Zhang   PetscErrorCode ierr;
941814e85d6SHong Zhang 
942814e85d6SHong Zhang   PetscFunctionBegin;
943814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
944c9ad9de0SHong Zhang   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
945814e85d6SHong Zhang 
946814e85d6SHong Zhang   PetscStackPush("TS user DRDY function for sensitivity analysis");
947c9ad9de0SHong Zhang   ierr = (*ts->drdufunction)(ts,t,U,DRDU,ts->costintegrandctx); CHKERRQ(ierr);
948814e85d6SHong Zhang   PetscStackPop;
949814e85d6SHong Zhang   PetscFunctionReturn(0);
950814e85d6SHong Zhang }
951814e85d6SHong Zhang 
952814e85d6SHong Zhang /*@
953814e85d6SHong Zhang   TSAdjointComputeDRDPFunction - Deprecated, use TSComputeDRDPFunction()
954814e85d6SHong Zhang 
955814e85d6SHong Zhang   Level: deprecated
956814e85d6SHong Zhang 
957814e85d6SHong Zhang @*/
958c9ad9de0SHong Zhang PetscErrorCode TSAdjointComputeDRDPFunction(TS ts,PetscReal t,Vec U,Vec *DRDP)
959814e85d6SHong Zhang {
960814e85d6SHong Zhang   PetscErrorCode ierr;
961814e85d6SHong Zhang 
962814e85d6SHong Zhang   PetscFunctionBegin;
963814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
964c9ad9de0SHong Zhang   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
965814e85d6SHong Zhang 
966814e85d6SHong Zhang   PetscStackPush("TS user DRDP function for sensitivity analysis");
967c9ad9de0SHong Zhang   ierr = (*ts->drdpfunction)(ts,t,U,DRDP,ts->costintegrandctx); CHKERRQ(ierr);
968814e85d6SHong Zhang   PetscStackPop;
969814e85d6SHong Zhang   PetscFunctionReturn(0);
970814e85d6SHong Zhang }
971814e85d6SHong Zhang 
972814e85d6SHong Zhang /*@C
973814e85d6SHong Zhang    TSAdjointMonitorSensi - monitors the first lambda sensitivity
974814e85d6SHong Zhang 
975814e85d6SHong Zhang    Level: intermediate
976814e85d6SHong Zhang 
977814e85d6SHong Zhang .keywords: TS, set, monitor
978814e85d6SHong Zhang 
979814e85d6SHong Zhang .seealso: TSAdjointMonitorSet()
980814e85d6SHong Zhang @*/
981814e85d6SHong Zhang PetscErrorCode TSAdjointMonitorSensi(TS ts,PetscInt step,PetscReal ptime,Vec v,PetscInt numcost,Vec *lambda,Vec *mu,PetscViewerAndFormat *vf)
982814e85d6SHong Zhang {
983814e85d6SHong Zhang   PetscErrorCode ierr;
984814e85d6SHong Zhang   PetscViewer    viewer = vf->viewer;
985814e85d6SHong Zhang 
986814e85d6SHong Zhang   PetscFunctionBegin;
987814e85d6SHong Zhang   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4);
988814e85d6SHong Zhang   ierr = PetscViewerPushFormat(viewer,vf->format);CHKERRQ(ierr);
989814e85d6SHong Zhang   ierr = VecView(lambda[0],viewer);CHKERRQ(ierr);
990814e85d6SHong Zhang   ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
991814e85d6SHong Zhang   PetscFunctionReturn(0);
992814e85d6SHong Zhang }
993814e85d6SHong Zhang 
994814e85d6SHong Zhang /*@C
995814e85d6SHong Zhang    TSAdjointMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user
996814e85d6SHong Zhang 
997814e85d6SHong Zhang    Collective on TS
998814e85d6SHong Zhang 
999814e85d6SHong Zhang    Input Parameters:
1000814e85d6SHong Zhang +  ts - TS object you wish to monitor
1001814e85d6SHong Zhang .  name - the monitor type one is seeking
1002814e85d6SHong Zhang .  help - message indicating what monitoring is done
1003814e85d6SHong Zhang .  manual - manual page for the monitor
1004814e85d6SHong Zhang .  monitor - the monitor function
1005814e85d6SHong Zhang -  monitorsetup - a function that is called once ONLY if the user selected this monitor that may set additional features of the TS or PetscViewer objects
1006814e85d6SHong Zhang 
1007814e85d6SHong Zhang    Level: developer
1008814e85d6SHong Zhang 
1009814e85d6SHong Zhang .seealso: PetscOptionsGetViewer(), PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(),
1010814e85d6SHong Zhang           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
1011814e85d6SHong Zhang           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
1012814e85d6SHong Zhang           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
1013814e85d6SHong Zhang           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1014814e85d6SHong Zhang           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1015814e85d6SHong Zhang           PetscOptionsFList(), PetscOptionsEList()
1016814e85d6SHong Zhang @*/
1017814e85d6SHong Zhang PetscErrorCode TSAdjointMonitorSetFromOptions(TS ts,const char name[],const char help[], const char manual[],PetscErrorCode (*monitor)(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,PetscViewerAndFormat*),PetscErrorCode (*monitorsetup)(TS,PetscViewerAndFormat*))
1018814e85d6SHong Zhang {
1019814e85d6SHong Zhang   PetscErrorCode    ierr;
1020814e85d6SHong Zhang   PetscViewer       viewer;
1021814e85d6SHong Zhang   PetscViewerFormat format;
1022814e85d6SHong Zhang   PetscBool         flg;
1023814e85d6SHong Zhang 
1024814e85d6SHong Zhang   PetscFunctionBegin;
102516413a6aSBarry Smith   ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)ts),((PetscObject) ts)->options,((PetscObject)ts)->prefix,name,&viewer,&format,&flg);CHKERRQ(ierr);
1026814e85d6SHong Zhang   if (flg) {
1027814e85d6SHong Zhang     PetscViewerAndFormat *vf;
1028814e85d6SHong Zhang     ierr = PetscViewerAndFormatCreate(viewer,format,&vf);CHKERRQ(ierr);
1029814e85d6SHong Zhang     ierr = PetscObjectDereference((PetscObject)viewer);CHKERRQ(ierr);
1030814e85d6SHong Zhang     if (monitorsetup) {
1031814e85d6SHong Zhang       ierr = (*monitorsetup)(ts,vf);CHKERRQ(ierr);
1032814e85d6SHong Zhang     }
1033814e85d6SHong Zhang     ierr = TSAdjointMonitorSet(ts,(PetscErrorCode (*)(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*))monitor,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy);CHKERRQ(ierr);
1034814e85d6SHong Zhang   }
1035814e85d6SHong Zhang   PetscFunctionReturn(0);
1036814e85d6SHong Zhang }
1037814e85d6SHong Zhang 
1038814e85d6SHong Zhang /*@C
1039814e85d6SHong Zhang    TSAdjointMonitorSet - Sets an ADDITIONAL function that is to be used at every
1040814e85d6SHong Zhang    timestep to display the iteration's  progress.
1041814e85d6SHong Zhang 
1042814e85d6SHong Zhang    Logically Collective on TS
1043814e85d6SHong Zhang 
1044814e85d6SHong Zhang    Input Parameters:
1045814e85d6SHong Zhang +  ts - the TS context obtained from TSCreate()
1046814e85d6SHong Zhang .  adjointmonitor - monitoring routine
1047814e85d6SHong Zhang .  adjointmctx - [optional] user-defined context for private data for the
1048814e85d6SHong Zhang              monitor routine (use NULL if no context is desired)
1049814e85d6SHong Zhang -  adjointmonitordestroy - [optional] routine that frees monitor context
1050814e85d6SHong Zhang           (may be NULL)
1051814e85d6SHong Zhang 
1052814e85d6SHong Zhang    Calling sequence of monitor:
1053814e85d6SHong Zhang $    int adjointmonitor(TS ts,PetscInt steps,PetscReal time,Vec u,PetscInt numcost,Vec *lambda, Vec *mu,void *adjointmctx)
1054814e85d6SHong Zhang 
1055814e85d6SHong Zhang +    ts - the TS context
1056814e85d6SHong Zhang .    steps - iteration number (after the final time step the monitor routine is called with a step of -1, this is at the final time which may have
1057814e85d6SHong Zhang                                been interpolated to)
1058814e85d6SHong Zhang .    time - current time
1059814e85d6SHong Zhang .    u - current iterate
1060814e85d6SHong Zhang .    numcost - number of cost functionos
1061814e85d6SHong Zhang .    lambda - sensitivities to initial conditions
1062814e85d6SHong Zhang .    mu - sensitivities to parameters
1063814e85d6SHong Zhang -    adjointmctx - [optional] adjoint monitoring context
1064814e85d6SHong Zhang 
1065814e85d6SHong Zhang    Notes:
1066814e85d6SHong Zhang    This routine adds an additional monitor to the list of monitors that
1067814e85d6SHong Zhang    already has been loaded.
1068814e85d6SHong Zhang 
1069814e85d6SHong Zhang    Fortran Notes:
1070814e85d6SHong Zhang     Only a single monitor function can be set for each TS object
1071814e85d6SHong Zhang 
1072814e85d6SHong Zhang    Level: intermediate
1073814e85d6SHong Zhang 
1074814e85d6SHong Zhang .keywords: TS, timestep, set, adjoint, monitor
1075814e85d6SHong Zhang 
1076814e85d6SHong Zhang .seealso: TSAdjointMonitorCancel()
1077814e85d6SHong Zhang @*/
1078814e85d6SHong Zhang PetscErrorCode TSAdjointMonitorSet(TS ts,PetscErrorCode (*adjointmonitor)(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*),void *adjointmctx,PetscErrorCode (*adjointmdestroy)(void**))
1079814e85d6SHong Zhang {
1080814e85d6SHong Zhang   PetscErrorCode ierr;
1081814e85d6SHong Zhang   PetscInt       i;
1082814e85d6SHong Zhang   PetscBool      identical;
1083814e85d6SHong Zhang 
1084814e85d6SHong Zhang   PetscFunctionBegin;
1085814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1086814e85d6SHong Zhang   for (i=0; i<ts->numbermonitors;i++) {
1087814e85d6SHong Zhang     ierr = PetscMonitorCompare((PetscErrorCode (*)(void))adjointmonitor,adjointmctx,adjointmdestroy,(PetscErrorCode (*)(void))ts->adjointmonitor[i],ts->adjointmonitorcontext[i],ts->adjointmonitordestroy[i],&identical);CHKERRQ(ierr);
1088814e85d6SHong Zhang     if (identical) PetscFunctionReturn(0);
1089814e85d6SHong Zhang   }
1090814e85d6SHong Zhang   if (ts->numberadjointmonitors >= MAXTSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many adjoint monitors set");
1091814e85d6SHong Zhang   ts->adjointmonitor[ts->numberadjointmonitors]          = adjointmonitor;
1092814e85d6SHong Zhang   ts->adjointmonitordestroy[ts->numberadjointmonitors]   = adjointmdestroy;
1093814e85d6SHong Zhang   ts->adjointmonitorcontext[ts->numberadjointmonitors++] = (void*)adjointmctx;
1094814e85d6SHong Zhang   PetscFunctionReturn(0);
1095814e85d6SHong Zhang }
1096814e85d6SHong Zhang 
1097814e85d6SHong Zhang /*@C
1098814e85d6SHong Zhang    TSAdjointMonitorCancel - Clears all the adjoint monitors that have been set on a time-step object.
1099814e85d6SHong Zhang 
1100814e85d6SHong Zhang    Logically Collective on TS
1101814e85d6SHong Zhang 
1102814e85d6SHong Zhang    Input Parameters:
1103814e85d6SHong Zhang .  ts - the TS context obtained from TSCreate()
1104814e85d6SHong Zhang 
1105814e85d6SHong Zhang    Notes:
1106814e85d6SHong Zhang    There is no way to remove a single, specific monitor.
1107814e85d6SHong Zhang 
1108814e85d6SHong Zhang    Level: intermediate
1109814e85d6SHong Zhang 
1110814e85d6SHong Zhang .keywords: TS, timestep, set, adjoint, monitor
1111814e85d6SHong Zhang 
1112814e85d6SHong Zhang .seealso: TSAdjointMonitorSet()
1113814e85d6SHong Zhang @*/
1114814e85d6SHong Zhang PetscErrorCode TSAdjointMonitorCancel(TS ts)
1115814e85d6SHong Zhang {
1116814e85d6SHong Zhang   PetscErrorCode ierr;
1117814e85d6SHong Zhang   PetscInt       i;
1118814e85d6SHong Zhang 
1119814e85d6SHong Zhang   PetscFunctionBegin;
1120814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1121814e85d6SHong Zhang   for (i=0; i<ts->numberadjointmonitors; i++) {
1122814e85d6SHong Zhang     if (ts->adjointmonitordestroy[i]) {
1123814e85d6SHong Zhang       ierr = (*ts->adjointmonitordestroy[i])(&ts->adjointmonitorcontext[i]);CHKERRQ(ierr);
1124814e85d6SHong Zhang     }
1125814e85d6SHong Zhang   }
1126814e85d6SHong Zhang   ts->numberadjointmonitors = 0;
1127814e85d6SHong Zhang   PetscFunctionReturn(0);
1128814e85d6SHong Zhang }
1129814e85d6SHong Zhang 
1130814e85d6SHong Zhang /*@C
1131814e85d6SHong Zhang    TSAdjointMonitorDefault - the default monitor of adjoint computations
1132814e85d6SHong Zhang 
1133814e85d6SHong Zhang    Level: intermediate
1134814e85d6SHong Zhang 
1135814e85d6SHong Zhang .keywords: TS, set, monitor
1136814e85d6SHong Zhang 
1137814e85d6SHong Zhang .seealso: TSAdjointMonitorSet()
1138814e85d6SHong Zhang @*/
1139814e85d6SHong Zhang PetscErrorCode TSAdjointMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,PetscInt numcost,Vec *lambda,Vec *mu,PetscViewerAndFormat *vf)
1140814e85d6SHong Zhang {
1141814e85d6SHong Zhang   PetscErrorCode ierr;
1142814e85d6SHong Zhang   PetscViewer    viewer = vf->viewer;
1143814e85d6SHong Zhang 
1144814e85d6SHong Zhang   PetscFunctionBegin;
1145814e85d6SHong Zhang   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4);
1146814e85d6SHong Zhang   ierr = PetscViewerPushFormat(viewer,vf->format);CHKERRQ(ierr);
1147814e85d6SHong Zhang   ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
1148814e85d6SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"%D TS dt %g time %g%s",step,(double)ts->time_step,(double)ptime,ts->steprollback ? " (r)\n" : "\n");CHKERRQ(ierr);
1149814e85d6SHong Zhang   ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
1150814e85d6SHong Zhang   ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
1151814e85d6SHong Zhang   PetscFunctionReturn(0);
1152814e85d6SHong Zhang }
1153814e85d6SHong Zhang 
1154814e85d6SHong Zhang /*@C
1155814e85d6SHong Zhang    TSAdjointMonitorDrawSensi - Monitors progress of the adjoint TS solvers by calling
1156814e85d6SHong Zhang    VecView() for the sensitivities to initial states at each timestep
1157814e85d6SHong Zhang 
1158814e85d6SHong Zhang    Collective on TS
1159814e85d6SHong Zhang 
1160814e85d6SHong Zhang    Input Parameters:
1161814e85d6SHong Zhang +  ts - the TS context
1162814e85d6SHong Zhang .  step - current time-step
1163814e85d6SHong Zhang .  ptime - current time
1164814e85d6SHong Zhang .  u - current state
1165814e85d6SHong Zhang .  numcost - number of cost functions
1166814e85d6SHong Zhang .  lambda - sensitivities to initial conditions
1167814e85d6SHong Zhang .  mu - sensitivities to parameters
1168814e85d6SHong Zhang -  dummy - either a viewer or NULL
1169814e85d6SHong Zhang 
1170814e85d6SHong Zhang    Level: intermediate
1171814e85d6SHong Zhang 
1172814e85d6SHong Zhang .keywords: TS,  vector, adjoint, monitor, view
1173814e85d6SHong Zhang 
1174814e85d6SHong Zhang .seealso: TSAdjointMonitorSet(), TSAdjointMonitorDefault(), VecView()
1175814e85d6SHong Zhang @*/
1176814e85d6SHong Zhang PetscErrorCode TSAdjointMonitorDrawSensi(TS ts,PetscInt step,PetscReal ptime,Vec u,PetscInt numcost,Vec *lambda,Vec *mu,void *dummy)
1177814e85d6SHong Zhang {
1178814e85d6SHong Zhang   PetscErrorCode   ierr;
1179814e85d6SHong Zhang   TSMonitorDrawCtx ictx = (TSMonitorDrawCtx)dummy;
1180814e85d6SHong Zhang   PetscDraw        draw;
1181814e85d6SHong Zhang   PetscReal        xl,yl,xr,yr,h;
1182814e85d6SHong Zhang   char             time[32];
1183814e85d6SHong Zhang 
1184814e85d6SHong Zhang   PetscFunctionBegin;
1185814e85d6SHong Zhang   if (!(((ictx->howoften > 0) && (!(step % ictx->howoften))) || ((ictx->howoften == -1) && ts->reason))) PetscFunctionReturn(0);
1186814e85d6SHong Zhang 
1187814e85d6SHong Zhang   ierr = VecView(lambda[0],ictx->viewer);CHKERRQ(ierr);
1188814e85d6SHong Zhang   ierr = PetscViewerDrawGetDraw(ictx->viewer,0,&draw);CHKERRQ(ierr);
1189814e85d6SHong Zhang   ierr = PetscSNPrintf(time,32,"Timestep %d Time %g",(int)step,(double)ptime);CHKERRQ(ierr);
1190814e85d6SHong Zhang   ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
1191814e85d6SHong Zhang   h    = yl + .95*(yr - yl);
1192814e85d6SHong Zhang   ierr = PetscDrawStringCentered(draw,.5*(xl+xr),h,PETSC_DRAW_BLACK,time);CHKERRQ(ierr);
1193814e85d6SHong Zhang   ierr = PetscDrawFlush(draw);CHKERRQ(ierr);
1194814e85d6SHong Zhang   PetscFunctionReturn(0);
1195814e85d6SHong Zhang }
1196814e85d6SHong Zhang 
1197814e85d6SHong Zhang /*
1198814e85d6SHong Zhang    TSAdjointSetFromOptions - Sets various TSAdjoint parameters from user options.
1199814e85d6SHong Zhang 
1200814e85d6SHong Zhang    Collective on TSAdjoint
1201814e85d6SHong Zhang 
1202814e85d6SHong Zhang    Input Parameter:
1203814e85d6SHong Zhang .  ts - the TS context
1204814e85d6SHong Zhang 
1205814e85d6SHong Zhang    Options Database Keys:
1206814e85d6SHong Zhang +  -ts_adjoint_solve <yes,no> After solving the ODE/DAE solve the adjoint problem (requires -ts_save_trajectory)
1207814e85d6SHong Zhang .  -ts_adjoint_monitor - print information at each adjoint time step
1208814e85d6SHong Zhang -  -ts_adjoint_monitor_draw_sensi - monitor the sensitivity of the first cost function wrt initial conditions (lambda[0]) graphically
1209814e85d6SHong Zhang 
1210814e85d6SHong Zhang    Level: developer
1211814e85d6SHong Zhang 
1212814e85d6SHong Zhang    Notes:
1213814e85d6SHong Zhang     This is not normally called directly by users
1214814e85d6SHong Zhang 
1215814e85d6SHong Zhang .keywords: TS, trajectory, timestep, set, options, database
1216814e85d6SHong Zhang 
1217814e85d6SHong Zhang .seealso: TSSetSaveTrajectory(), TSTrajectorySetUp()
1218814e85d6SHong Zhang */
1219814e85d6SHong Zhang PetscErrorCode TSAdjointSetFromOptions(PetscOptionItems *PetscOptionsObject,TS ts)
1220814e85d6SHong Zhang {
1221814e85d6SHong Zhang   PetscBool      tflg,opt;
1222814e85d6SHong Zhang   PetscErrorCode ierr;
1223814e85d6SHong Zhang 
1224814e85d6SHong Zhang   PetscFunctionBegin;
1225814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,2);
1226814e85d6SHong Zhang   ierr = PetscOptionsHead(PetscOptionsObject,"TS Adjoint options");CHKERRQ(ierr);
1227814e85d6SHong Zhang   tflg = ts->adjoint_solve ? PETSC_TRUE : PETSC_FALSE;
1228814e85d6SHong Zhang   ierr = PetscOptionsBool("-ts_adjoint_solve","Solve the adjoint problem immediately after solving the forward problem","",tflg,&tflg,&opt);CHKERRQ(ierr);
1229814e85d6SHong Zhang   if (opt) {
1230814e85d6SHong Zhang     ierr = TSSetSaveTrajectory(ts);CHKERRQ(ierr);
1231814e85d6SHong Zhang     ts->adjoint_solve = tflg;
1232814e85d6SHong Zhang   }
1233814e85d6SHong Zhang   ierr = TSAdjointMonitorSetFromOptions(ts,"-ts_adjoint_monitor","Monitor adjoint timestep size","TSAdjointMonitorDefault",TSAdjointMonitorDefault,NULL);CHKERRQ(ierr);
1234814e85d6SHong Zhang   ierr = TSAdjointMonitorSetFromOptions(ts,"-ts_adjoint_monitor_sensi","Monitor sensitivity in the adjoint computation","TSAdjointMonitorSensi",TSAdjointMonitorSensi,NULL);CHKERRQ(ierr);
1235814e85d6SHong Zhang   opt  = PETSC_FALSE;
1236814e85d6SHong Zhang   ierr = PetscOptionsName("-ts_adjoint_monitor_draw_sensi","Monitor adjoint sensitivities (lambda only) graphically","TSAdjointMonitorDrawSensi",&opt);CHKERRQ(ierr);
1237814e85d6SHong Zhang   if (opt) {
1238814e85d6SHong Zhang     TSMonitorDrawCtx ctx;
1239814e85d6SHong Zhang     PetscInt         howoften = 1;
1240814e85d6SHong Zhang 
1241814e85d6SHong Zhang     ierr = PetscOptionsInt("-ts_adjoint_monitor_draw_sensi","Monitor adjoint sensitivities (lambda only) graphically","TSAdjointMonitorDrawSensi",howoften,&howoften,NULL);CHKERRQ(ierr);
1242814e85d6SHong Zhang     ierr = TSMonitorDrawCtxCreate(PetscObjectComm((PetscObject)ts),0,0,PETSC_DECIDE,PETSC_DECIDE,300,300,howoften,&ctx);CHKERRQ(ierr);
1243814e85d6SHong Zhang     ierr = TSAdjointMonitorSet(ts,TSAdjointMonitorDrawSensi,ctx,(PetscErrorCode (*)(void**))TSMonitorDrawCtxDestroy);CHKERRQ(ierr);
1244814e85d6SHong Zhang   }
1245814e85d6SHong Zhang   PetscFunctionReturn(0);
1246814e85d6SHong Zhang }
1247814e85d6SHong Zhang 
1248814e85d6SHong Zhang /*@
1249814e85d6SHong Zhang    TSAdjointStep - Steps one time step backward in the adjoint run
1250814e85d6SHong Zhang 
1251814e85d6SHong Zhang    Collective on TS
1252814e85d6SHong Zhang 
1253814e85d6SHong Zhang    Input Parameter:
1254814e85d6SHong Zhang .  ts - the TS context obtained from TSCreate()
1255814e85d6SHong Zhang 
1256814e85d6SHong Zhang    Level: intermediate
1257814e85d6SHong Zhang 
1258814e85d6SHong Zhang .keywords: TS, adjoint, step
1259814e85d6SHong Zhang 
1260814e85d6SHong Zhang .seealso: TSAdjointSetUp(), TSAdjointSolve()
1261814e85d6SHong Zhang @*/
1262814e85d6SHong Zhang PetscErrorCode TSAdjointStep(TS ts)
1263814e85d6SHong Zhang {
1264814e85d6SHong Zhang   DM               dm;
1265814e85d6SHong Zhang   PetscErrorCode   ierr;
1266814e85d6SHong Zhang 
1267814e85d6SHong Zhang   PetscFunctionBegin;
1268814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1269814e85d6SHong Zhang   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
1270814e85d6SHong Zhang   ierr = TSAdjointSetUp(ts);CHKERRQ(ierr);
1271814e85d6SHong Zhang 
1272814e85d6SHong Zhang   ts->reason = TS_CONVERGED_ITERATING;
1273814e85d6SHong Zhang   ts->ptime_prev = ts->ptime;
1274814e85d6SHong Zhang   if (!ts->ops->adjointstep) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_NOT_CONVERGED,"TSStep has failed because the adjoint of  %s has not been implemented, try other time stepping methods for adjoint sensitivity analysis",((PetscObject)ts)->type_name);
1275814e85d6SHong Zhang   ierr = PetscLogEventBegin(TS_AdjointStep,ts,0,0,0);CHKERRQ(ierr);
1276814e85d6SHong Zhang   ierr = (*ts->ops->adjointstep)(ts);CHKERRQ(ierr);
1277814e85d6SHong Zhang   ierr = PetscLogEventEnd(TS_AdjointStep,ts,0,0,0);CHKERRQ(ierr);
1278814e85d6SHong Zhang   ts->adjoint_steps++; ts->steps--;
1279814e85d6SHong Zhang 
1280814e85d6SHong Zhang   if (ts->reason < 0) {
1281814e85d6SHong Zhang     if (ts->errorifstepfailed) {
128205c9950eSHong Zhang       if (ts->reason == TSADJOINT_DIVERGED_LINEAR_SOLVE) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_NOT_CONVERGED,"TSAdjointStep has failed due to %s",TSConvergedReasons[ts->reason]);
128305c9950eSHong Zhang       else SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_NOT_CONVERGED,"TSAdjointStep has failed due to %s",TSConvergedReasons[ts->reason]);
1284814e85d6SHong Zhang     }
1285814e85d6SHong Zhang   } else if (!ts->reason) {
1286814e85d6SHong Zhang     if (ts->adjoint_steps >= ts->adjoint_max_steps) ts->reason = TS_CONVERGED_ITS;
1287814e85d6SHong Zhang   }
1288814e85d6SHong Zhang   PetscFunctionReturn(0);
1289814e85d6SHong Zhang }
1290814e85d6SHong Zhang 
1291814e85d6SHong Zhang /*@
1292814e85d6SHong Zhang    TSAdjointSolve - Solves the discrete ajoint problem for an ODE/DAE
1293814e85d6SHong Zhang 
1294814e85d6SHong Zhang    Collective on TS
1295814e85d6SHong Zhang 
1296814e85d6SHong Zhang    Input Parameter:
1297814e85d6SHong Zhang .  ts - the TS context obtained from TSCreate()
1298814e85d6SHong Zhang 
1299814e85d6SHong Zhang    Options Database:
1300814e85d6SHong Zhang . -ts_adjoint_view_solution <viewerinfo> - views the first gradient with respect to the initial values
1301814e85d6SHong Zhang 
1302814e85d6SHong Zhang    Level: intermediate
1303814e85d6SHong Zhang 
1304814e85d6SHong Zhang    Notes:
1305814e85d6SHong Zhang    This must be called after a call to TSSolve() that solves the forward problem
1306814e85d6SHong Zhang 
1307814e85d6SHong Zhang    By default this will integrate back to the initial time, one can use TSAdjointSetSteps() to step back to a later time
1308814e85d6SHong Zhang 
1309814e85d6SHong Zhang .keywords: TS, timestep, solve
1310814e85d6SHong Zhang 
1311814e85d6SHong Zhang .seealso: TSCreate(), TSSetCostGradients(), TSSetSolution(), TSAdjointStep()
1312814e85d6SHong Zhang @*/
1313814e85d6SHong Zhang PetscErrorCode TSAdjointSolve(TS ts)
1314814e85d6SHong Zhang {
1315814e85d6SHong Zhang   PetscErrorCode    ierr;
1316814e85d6SHong Zhang 
1317814e85d6SHong Zhang   PetscFunctionBegin;
1318814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1319814e85d6SHong Zhang   ierr = TSAdjointSetUp(ts);CHKERRQ(ierr);
1320814e85d6SHong Zhang 
1321814e85d6SHong Zhang   /* reset time step and iteration counters */
1322814e85d6SHong Zhang   ts->adjoint_steps     = 0;
1323814e85d6SHong Zhang   ts->ksp_its           = 0;
1324814e85d6SHong Zhang   ts->snes_its          = 0;
1325814e85d6SHong Zhang   ts->num_snes_failures = 0;
1326814e85d6SHong Zhang   ts->reject            = 0;
1327814e85d6SHong Zhang   ts->reason            = TS_CONVERGED_ITERATING;
1328814e85d6SHong Zhang 
1329814e85d6SHong Zhang   if (!ts->adjoint_max_steps) ts->adjoint_max_steps = ts->steps;
1330814e85d6SHong Zhang   if (ts->adjoint_steps >= ts->adjoint_max_steps) ts->reason = TS_CONVERGED_ITS;
1331814e85d6SHong Zhang 
1332814e85d6SHong Zhang   while (!ts->reason) {
1333814e85d6SHong Zhang     ierr = TSTrajectoryGet(ts->trajectory,ts,ts->steps,&ts->ptime);CHKERRQ(ierr);
1334814e85d6SHong Zhang     ierr = TSAdjointMonitor(ts,ts->steps,ts->ptime,ts->vec_sol,ts->numcost,ts->vecs_sensi,ts->vecs_sensip);CHKERRQ(ierr);
1335814e85d6SHong Zhang     ierr = TSAdjointEventHandler(ts);CHKERRQ(ierr);
1336814e85d6SHong Zhang     ierr = TSAdjointStep(ts);CHKERRQ(ierr);
1337814e85d6SHong Zhang     if (ts->vec_costintegral && !ts->costintegralfwd) {
1338814e85d6SHong Zhang       ierr = TSAdjointCostIntegral(ts);CHKERRQ(ierr);
1339814e85d6SHong Zhang     }
1340814e85d6SHong Zhang   }
1341814e85d6SHong Zhang   ierr = TSTrajectoryGet(ts->trajectory,ts,ts->steps,&ts->ptime);CHKERRQ(ierr);
1342814e85d6SHong Zhang   ierr = TSAdjointMonitor(ts,ts->steps,ts->ptime,ts->vec_sol,ts->numcost,ts->vecs_sensi,ts->vecs_sensip);CHKERRQ(ierr);
1343814e85d6SHong Zhang   ts->solvetime = ts->ptime;
1344814e85d6SHong Zhang   ierr = TSTrajectoryViewFromOptions(ts->trajectory,NULL,"-ts_trajectory_view");CHKERRQ(ierr);
1345814e85d6SHong Zhang   ierr = VecViewFromOptions(ts->vecs_sensi[0],(PetscObject) ts, "-ts_adjoint_view_solution");CHKERRQ(ierr);
1346814e85d6SHong Zhang   ts->adjoint_max_steps = 0;
1347814e85d6SHong Zhang   PetscFunctionReturn(0);
1348814e85d6SHong Zhang }
1349814e85d6SHong Zhang 
1350814e85d6SHong Zhang /*@C
1351814e85d6SHong Zhang    TSAdjointMonitor - Runs all user-provided adjoint monitor routines set using TSAdjointMonitorSet()
1352814e85d6SHong Zhang 
1353814e85d6SHong Zhang    Collective on TS
1354814e85d6SHong Zhang 
1355814e85d6SHong Zhang    Input Parameters:
1356814e85d6SHong Zhang +  ts - time stepping context obtained from TSCreate()
1357814e85d6SHong Zhang .  step - step number that has just completed
1358814e85d6SHong Zhang .  ptime - model time of the state
1359814e85d6SHong Zhang .  u - state at the current model time
1360814e85d6SHong Zhang .  numcost - number of cost functions (dimension of lambda  or mu)
1361814e85d6SHong Zhang .  lambda - vectors containing the gradients of the cost functions with respect to the ODE/DAE solution variables
1362814e85d6SHong Zhang -  mu - vectors containing the gradients of the cost functions with respect to the problem parameters
1363814e85d6SHong Zhang 
1364814e85d6SHong Zhang    Notes:
1365814e85d6SHong Zhang    TSAdjointMonitor() is typically used automatically within the time stepping implementations.
1366814e85d6SHong Zhang    Users would almost never call this routine directly.
1367814e85d6SHong Zhang 
1368814e85d6SHong Zhang    Level: developer
1369814e85d6SHong Zhang 
1370814e85d6SHong Zhang .keywords: TS, timestep
1371814e85d6SHong Zhang @*/
1372814e85d6SHong Zhang PetscErrorCode TSAdjointMonitor(TS ts,PetscInt step,PetscReal ptime,Vec u,PetscInt numcost,Vec *lambda, Vec *mu)
1373814e85d6SHong Zhang {
1374814e85d6SHong Zhang   PetscErrorCode ierr;
1375814e85d6SHong Zhang   PetscInt       i,n = ts->numberadjointmonitors;
1376814e85d6SHong Zhang 
1377814e85d6SHong Zhang   PetscFunctionBegin;
1378814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1379814e85d6SHong Zhang   PetscValidHeaderSpecific(u,VEC_CLASSID,4);
13808860a134SJunchao Zhang   ierr = VecLockReadPush(u);CHKERRQ(ierr);
1381814e85d6SHong Zhang   for (i=0; i<n; i++) {
1382814e85d6SHong Zhang     ierr = (*ts->adjointmonitor[i])(ts,step,ptime,u,numcost,lambda,mu,ts->adjointmonitorcontext[i]);CHKERRQ(ierr);
1383814e85d6SHong Zhang   }
13848860a134SJunchao Zhang   ierr = VecLockReadPop(u);CHKERRQ(ierr);
1385814e85d6SHong Zhang   PetscFunctionReturn(0);
1386814e85d6SHong Zhang }
1387814e85d6SHong Zhang 
1388814e85d6SHong Zhang /*@
1389814e85d6SHong Zhang  TSAdjointCostIntegral - Evaluate the cost integral in the adjoint run.
1390814e85d6SHong Zhang 
1391814e85d6SHong Zhang  Collective on TS
1392814e85d6SHong Zhang 
1393814e85d6SHong Zhang  Input Arguments:
1394814e85d6SHong Zhang  .  ts - time stepping context
1395814e85d6SHong Zhang 
1396814e85d6SHong Zhang  Level: advanced
1397814e85d6SHong Zhang 
1398814e85d6SHong Zhang  Notes:
1399814e85d6SHong Zhang  This function cannot be called until TSAdjointStep() has been completed.
1400814e85d6SHong Zhang 
1401814e85d6SHong Zhang  .seealso: TSAdjointSolve(), TSAdjointStep
1402814e85d6SHong Zhang  @*/
1403814e85d6SHong Zhang PetscErrorCode TSAdjointCostIntegral(TS ts)
1404814e85d6SHong Zhang {
1405814e85d6SHong Zhang     PetscErrorCode ierr;
1406814e85d6SHong Zhang     PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1407814e85d6SHong Zhang     if (!ts->ops->adjointintegral) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"%s does not provide integral evaluation in the adjoint run",((PetscObject)ts)->type_name);
1408814e85d6SHong Zhang     ierr = (*ts->ops->adjointintegral)(ts);CHKERRQ(ierr);
1409814e85d6SHong Zhang     PetscFunctionReturn(0);
1410814e85d6SHong Zhang }
1411814e85d6SHong Zhang 
1412814e85d6SHong Zhang /* ------------------ Forward (tangent linear) sensitivity  ------------------*/
1413814e85d6SHong Zhang 
1414814e85d6SHong Zhang /*@
1415814e85d6SHong Zhang   TSForwardSetUp - Sets up the internal data structures for the later use
1416814e85d6SHong Zhang   of forward sensitivity analysis
1417814e85d6SHong Zhang 
1418814e85d6SHong Zhang   Collective on TS
1419814e85d6SHong Zhang 
1420814e85d6SHong Zhang   Input Parameter:
1421814e85d6SHong Zhang . ts - the TS context obtained from TSCreate()
1422814e85d6SHong Zhang 
1423814e85d6SHong Zhang   Level: advanced
1424814e85d6SHong Zhang 
1425814e85d6SHong Zhang .keywords: TS, forward sensitivity, setup
1426814e85d6SHong Zhang 
1427814e85d6SHong Zhang .seealso: TSCreate(), TSDestroy(), TSSetUp()
1428814e85d6SHong Zhang @*/
1429814e85d6SHong Zhang PetscErrorCode TSForwardSetUp(TS ts)
1430814e85d6SHong Zhang {
1431814e85d6SHong Zhang   PetscErrorCode ierr;
1432814e85d6SHong Zhang 
1433814e85d6SHong Zhang   PetscFunctionBegin;
1434814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1435814e85d6SHong Zhang   if (ts->forwardsetupcalled) PetscFunctionReturn(0);
1436814e85d6SHong Zhang   if (ts->vec_costintegral && !ts->vecs_integral_sensip ) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TSForwardSetIntegralGradients() before TSSetCostIntegrand()");
1437814e85d6SHong Zhang   if (ts->vecs_integral_sensip) {
1438c9ad9de0SHong Zhang     ierr = VecDuplicateVecs(ts->vec_sol,ts->numcost,&ts->vecs_drdu);CHKERRQ(ierr);
1439814e85d6SHong Zhang     ierr = VecDuplicateVecs(ts->vecs_integral_sensip[0],ts->numcost,&ts->vecs_drdp);CHKERRQ(ierr);
1440814e85d6SHong Zhang   }
1441*33f72c5dSHong Zhang   if (!ts->Jacp && ts->Jacprhs) ts->Jacp = ts->Jacprhs;
1442814e85d6SHong Zhang 
1443814e85d6SHong Zhang   if (ts->ops->forwardsetup) {
1444814e85d6SHong Zhang     ierr = (*ts->ops->forwardsetup)(ts);CHKERRQ(ierr);
1445814e85d6SHong Zhang   }
1446814e85d6SHong Zhang   ts->forwardsetupcalled = PETSC_TRUE;
1447814e85d6SHong Zhang   PetscFunctionReturn(0);
1448814e85d6SHong Zhang }
1449814e85d6SHong Zhang 
1450814e85d6SHong Zhang /*@
14517adebddeSHong Zhang   TSForwardReset - Reset the internal data structures used by forward sensitivity analysis
14527adebddeSHong Zhang 
14537adebddeSHong Zhang   Collective on TS
14547adebddeSHong Zhang 
14557adebddeSHong Zhang   Input Parameter:
14567adebddeSHong Zhang . ts - the TS context obtained from TSCreate()
14577adebddeSHong Zhang 
14587adebddeSHong Zhang   Level: advanced
14597adebddeSHong Zhang 
14607adebddeSHong Zhang .keywords: TS, forward sensitivity, reset
14617adebddeSHong Zhang 
14627adebddeSHong Zhang .seealso: TSCreate(), TSDestroy(), TSForwardSetUp()
14637adebddeSHong Zhang @*/
14647adebddeSHong Zhang PetscErrorCode TSForwardReset(TS ts)
14657adebddeSHong Zhang {
14667adebddeSHong Zhang   PetscErrorCode ierr;
14677adebddeSHong Zhang 
14687adebddeSHong Zhang   PetscFunctionBegin;
14697adebddeSHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
14707adebddeSHong Zhang   if (ts->ops->forwardreset) {
14717adebddeSHong Zhang     ierr = (*ts->ops->forwardreset)(ts);CHKERRQ(ierr);
14727adebddeSHong Zhang   }
14737adebddeSHong Zhang   ierr = MatDestroy(&ts->mat_sensip);CHKERRQ(ierr);
14747adebddeSHong Zhang   ts->vecs_integral_sensip = NULL;
14757adebddeSHong Zhang   ts->forward_solve        = PETSC_FALSE;
14767adebddeSHong Zhang   ts->forwardsetupcalled   = PETSC_FALSE;
14777adebddeSHong Zhang   PetscFunctionReturn(0);
14787adebddeSHong Zhang }
14797adebddeSHong Zhang 
14807adebddeSHong Zhang /*@
1481814e85d6SHong Zhang   TSForwardSetIntegralGradients - Set the vectors holding forward sensitivities of the integral term.
1482814e85d6SHong Zhang 
1483814e85d6SHong Zhang   Input Parameter:
1484814e85d6SHong Zhang . ts- the TS context obtained from TSCreate()
1485814e85d6SHong Zhang . numfwdint- number of integrals
1486814e85d6SHong Zhang . vp = the vectors containing the gradients for each integral w.r.t. parameters
1487814e85d6SHong Zhang 
1488814e85d6SHong Zhang   Level: intermediate
1489814e85d6SHong Zhang 
1490814e85d6SHong Zhang .keywords: TS, forward sensitivity
1491814e85d6SHong Zhang 
1492814e85d6SHong Zhang .seealso: TSForwardGetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep()
1493814e85d6SHong Zhang @*/
1494814e85d6SHong Zhang PetscErrorCode TSForwardSetIntegralGradients(TS ts,PetscInt numfwdint,Vec *vp)
1495814e85d6SHong Zhang {
1496814e85d6SHong Zhang   PetscFunctionBegin;
1497814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1498814e85d6SHong Zhang   if (ts->numcost && ts->numcost!=numfwdint) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"The number of cost functions (2rd parameter of TSSetCostIntegrand()) is inconsistent with the one set by TSSetCostIntegrand()");
1499814e85d6SHong Zhang   if (!ts->numcost) ts->numcost = numfwdint;
1500814e85d6SHong Zhang 
1501814e85d6SHong Zhang   ts->vecs_integral_sensip = vp;
1502814e85d6SHong Zhang   PetscFunctionReturn(0);
1503814e85d6SHong Zhang }
1504814e85d6SHong Zhang 
1505814e85d6SHong Zhang /*@
1506814e85d6SHong Zhang   TSForwardGetIntegralGradients - Returns the forward sensitivities ofthe integral term.
1507814e85d6SHong Zhang 
1508814e85d6SHong Zhang   Input Parameter:
1509814e85d6SHong Zhang . ts- the TS context obtained from TSCreate()
1510814e85d6SHong Zhang 
1511814e85d6SHong Zhang   Output Parameter:
1512814e85d6SHong Zhang . vp = the vectors containing the gradients for each integral w.r.t. parameters
1513814e85d6SHong Zhang 
1514814e85d6SHong Zhang   Level: intermediate
1515814e85d6SHong Zhang 
1516814e85d6SHong Zhang .keywords: TS, forward sensitivity
1517814e85d6SHong Zhang 
1518814e85d6SHong Zhang .seealso: TSForwardSetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep()
1519814e85d6SHong Zhang @*/
1520814e85d6SHong Zhang PetscErrorCode TSForwardGetIntegralGradients(TS ts,PetscInt *numfwdint,Vec **vp)
1521814e85d6SHong Zhang {
1522814e85d6SHong Zhang   PetscFunctionBegin;
1523814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1524814e85d6SHong Zhang   PetscValidPointer(vp,3);
1525814e85d6SHong Zhang   if (numfwdint) *numfwdint = ts->numcost;
1526814e85d6SHong Zhang   if (vp) *vp = ts->vecs_integral_sensip;
1527814e85d6SHong Zhang   PetscFunctionReturn(0);
1528814e85d6SHong Zhang }
1529814e85d6SHong Zhang 
1530814e85d6SHong Zhang /*@
1531814e85d6SHong Zhang   TSForwardStep - Compute the forward sensitivity for one time step.
1532814e85d6SHong Zhang 
1533814e85d6SHong Zhang   Collective on TS
1534814e85d6SHong Zhang 
1535814e85d6SHong Zhang   Input Arguments:
1536814e85d6SHong Zhang . ts - time stepping context
1537814e85d6SHong Zhang 
1538814e85d6SHong Zhang   Level: advanced
1539814e85d6SHong Zhang 
1540814e85d6SHong Zhang   Notes:
1541814e85d6SHong Zhang   This function cannot be called until TSStep() has been completed.
1542814e85d6SHong Zhang 
1543814e85d6SHong Zhang .keywords: TS, forward sensitivity
1544814e85d6SHong Zhang 
1545814e85d6SHong Zhang .seealso: TSForwardSetSensitivities(), TSForwardGetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardSetUp()
1546814e85d6SHong Zhang @*/
1547814e85d6SHong Zhang PetscErrorCode TSForwardStep(TS ts)
1548814e85d6SHong Zhang {
1549814e85d6SHong Zhang   PetscErrorCode ierr;
1550814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1551814e85d6SHong Zhang   if (!ts->ops->forwardstep) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"%s does not provide forward sensitivity analysis",((PetscObject)ts)->type_name);
1552814e85d6SHong Zhang   ierr = PetscLogEventBegin(TS_ForwardStep,ts,0,0,0);CHKERRQ(ierr);
1553814e85d6SHong Zhang   ierr = (*ts->ops->forwardstep)(ts);CHKERRQ(ierr);
1554814e85d6SHong Zhang   ierr = PetscLogEventEnd(TS_ForwardStep,ts,0,0,0);CHKERRQ(ierr);
1555814e85d6SHong Zhang   PetscFunctionReturn(0);
1556814e85d6SHong Zhang }
1557814e85d6SHong Zhang 
1558814e85d6SHong Zhang /*@
1559814e85d6SHong Zhang   TSForwardSetSensitivities - Sets the initial value of the trajectory sensitivities of solution  w.r.t. the problem parameters and initial values.
1560814e85d6SHong Zhang 
1561814e85d6SHong Zhang   Logically Collective on TS and Vec
1562814e85d6SHong Zhang 
1563814e85d6SHong Zhang   Input Parameters:
1564814e85d6SHong Zhang + ts - the TS context obtained from TSCreate()
1565814e85d6SHong Zhang . nump - number of parameters
1566814e85d6SHong Zhang - Smat - sensitivities with respect to the parameters, the number of entries in these vectors is the same as the number of parameters
1567814e85d6SHong Zhang 
1568814e85d6SHong Zhang   Level: beginner
1569814e85d6SHong Zhang 
1570814e85d6SHong Zhang   Notes:
1571814e85d6SHong Zhang   Forward sensitivity is also called 'trajectory sensitivity' in some fields such as power systems.
1572814e85d6SHong Zhang   This function turns on a flag to trigger TSSolve() to compute forward sensitivities automatically.
1573814e85d6SHong Zhang   You must call this function before TSSolve().
1574814e85d6SHong Zhang   The entries in the sensitivity matrix must be correctly initialized with the values S = dy/dp|startingtime.
1575814e85d6SHong Zhang 
1576814e85d6SHong Zhang .keywords: TS, timestep, set, forward sensitivity, initial values
1577814e85d6SHong Zhang 
1578814e85d6SHong Zhang .seealso: TSForwardGetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep()
1579814e85d6SHong Zhang @*/
1580814e85d6SHong Zhang PetscErrorCode TSForwardSetSensitivities(TS ts,PetscInt nump,Mat Smat)
1581814e85d6SHong Zhang {
1582814e85d6SHong Zhang   PetscErrorCode ierr;
1583814e85d6SHong Zhang 
1584814e85d6SHong Zhang   PetscFunctionBegin;
1585814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1586814e85d6SHong Zhang   PetscValidHeaderSpecific(Smat,MAT_CLASSID,3);
1587814e85d6SHong Zhang   ts->forward_solve  = PETSC_TRUE;
1588b5b4017aSHong Zhang   if (nump == PETSC_DEFAULT) {
1589b5b4017aSHong Zhang     ierr = MatGetSize(Smat,NULL,&ts->num_parameters);CHKERRQ(ierr);
1590b5b4017aSHong Zhang   } else ts->num_parameters = nump;
1591814e85d6SHong Zhang   ierr = PetscObjectReference((PetscObject)Smat);CHKERRQ(ierr);
1592814e85d6SHong Zhang   ierr = MatDestroy(&ts->mat_sensip);CHKERRQ(ierr);
1593814e85d6SHong Zhang   ts->mat_sensip = Smat;
1594814e85d6SHong Zhang   PetscFunctionReturn(0);
1595814e85d6SHong Zhang }
1596814e85d6SHong Zhang 
1597814e85d6SHong Zhang /*@
1598814e85d6SHong Zhang   TSForwardGetSensitivities - Returns the trajectory sensitivities
1599814e85d6SHong Zhang 
1600814e85d6SHong Zhang   Not Collective, but Vec returned is parallel if TS is parallel
1601814e85d6SHong Zhang 
1602814e85d6SHong Zhang   Output Parameter:
1603814e85d6SHong Zhang + ts - the TS context obtained from TSCreate()
1604814e85d6SHong Zhang . nump - number of parameters
1605814e85d6SHong Zhang - Smat - sensitivities with respect to the parameters, the number of entries in these vectors is the same as the number of parameters
1606814e85d6SHong Zhang 
1607814e85d6SHong Zhang   Level: intermediate
1608814e85d6SHong Zhang 
1609814e85d6SHong Zhang .keywords: TS, forward sensitivity
1610814e85d6SHong Zhang 
1611814e85d6SHong Zhang .seealso: TSForwardSetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep()
1612814e85d6SHong Zhang @*/
1613814e85d6SHong Zhang PetscErrorCode TSForwardGetSensitivities(TS ts,PetscInt *nump,Mat *Smat)
1614814e85d6SHong Zhang {
1615814e85d6SHong Zhang   PetscFunctionBegin;
1616814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1617814e85d6SHong Zhang   if (nump) *nump = ts->num_parameters;
1618814e85d6SHong Zhang   if (Smat) *Smat = ts->mat_sensip;
1619814e85d6SHong Zhang   PetscFunctionReturn(0);
1620814e85d6SHong Zhang }
1621814e85d6SHong Zhang 
1622814e85d6SHong Zhang /*@
1623814e85d6SHong Zhang    TSForwardCostIntegral - Evaluate the cost integral in the forward run.
1624814e85d6SHong Zhang 
1625814e85d6SHong Zhang    Collective on TS
1626814e85d6SHong Zhang 
1627814e85d6SHong Zhang    Input Arguments:
1628814e85d6SHong Zhang .  ts - time stepping context
1629814e85d6SHong Zhang 
1630814e85d6SHong Zhang    Level: advanced
1631814e85d6SHong Zhang 
1632814e85d6SHong Zhang    Notes:
1633814e85d6SHong Zhang    This function cannot be called until TSStep() has been completed.
1634814e85d6SHong Zhang 
1635814e85d6SHong Zhang .seealso: TSSolve(), TSAdjointCostIntegral()
1636814e85d6SHong Zhang @*/
1637814e85d6SHong Zhang PetscErrorCode TSForwardCostIntegral(TS ts)
1638814e85d6SHong Zhang {
1639814e85d6SHong Zhang   PetscErrorCode ierr;
1640814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1641814e85d6SHong Zhang   if (!ts->ops->forwardintegral) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"%s does not provide integral evaluation in the forward run",((PetscObject)ts)->type_name);
1642814e85d6SHong Zhang   ierr = (*ts->ops->forwardintegral)(ts);CHKERRQ(ierr);
1643814e85d6SHong Zhang   PetscFunctionReturn(0);
1644814e85d6SHong Zhang }
1645b5b4017aSHong Zhang 
1646b5b4017aSHong Zhang /*@
1647b5b4017aSHong Zhang   TSForwardSetInitialSensitivities - Set initial values for tangent linear sensitivities
1648b5b4017aSHong Zhang 
1649b5b4017aSHong Zhang   Collective on TS and Mat
1650b5b4017aSHong Zhang 
1651b5b4017aSHong Zhang   Input Parameter
1652b5b4017aSHong Zhang + ts - the TS context obtained from TSCreate()
1653b5b4017aSHong Zhang - didp - parametric sensitivities of the initial condition
1654b5b4017aSHong Zhang 
1655b5b4017aSHong Zhang   Level: intermediate
1656b5b4017aSHong Zhang 
1657b5b4017aSHong Zhang   Notes: TSSolve() allows users to pass the initial solution directly to TS. But the tangent linear variables cannot be initialized in this way. This function is used to set initial values for tangent linear variables.
1658b5b4017aSHong Zhang 
1659b5b4017aSHong Zhang .seealso: TSForwardSetSensitivities()
1660b5b4017aSHong Zhang @*/
1661b5b4017aSHong Zhang PetscErrorCode TSForwardSetInitialSensitivities(TS ts,Mat didp)
1662b5b4017aSHong Zhang {
1663b5b4017aSHong Zhang   PetscErrorCode ierr;
1664b5b4017aSHong Zhang 
1665b5b4017aSHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1666b5b4017aSHong Zhang   PetscValidHeaderSpecific(didp,MAT_CLASSID,2);
1667b5b4017aSHong Zhang   if (!ts->mat_sensip) {
1668b5b4017aSHong Zhang     ierr = TSForwardSetSensitivities(ts,PETSC_DEFAULT,didp);CHKERRQ(ierr);
1669b5b4017aSHong Zhang   }
1670b5b4017aSHong Zhang   PetscFunctionReturn(0);
1671b5b4017aSHong Zhang }
16726affb6f8SHong Zhang 
16736affb6f8SHong Zhang /*@
16746affb6f8SHong Zhang    TSForwardGetStages - Get the number of stages and the tangent linear sensitivities at the intermediate stages
16756affb6f8SHong Zhang 
16766affb6f8SHong Zhang    Input Parameter:
16776affb6f8SHong Zhang .  ts - the TS context obtained from TSCreate()
16786affb6f8SHong Zhang 
16796affb6f8SHong Zhang    Output Parameters:
16806affb6f8SHong Zhang +  ns - nu
16816affb6f8SHong Zhang -  S - tangent linear sensitivities at the intermediate stages
16826affb6f8SHong Zhang 
16836affb6f8SHong Zhang    Level: advanced
16846affb6f8SHong Zhang 
16856affb6f8SHong Zhang .keywords: TS, second-order adjoint, forward sensitivity
16866affb6f8SHong Zhang @*/
16876affb6f8SHong Zhang PetscErrorCode TSForwardGetStages(TS ts,PetscInt *ns,Mat **S)
16886affb6f8SHong Zhang {
16896affb6f8SHong Zhang   PetscErrorCode ierr;
16906affb6f8SHong Zhang 
16916affb6f8SHong Zhang   PetscFunctionBegin;
16926affb6f8SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
16936affb6f8SHong Zhang 
16946affb6f8SHong Zhang   if (!ts->ops->getstages) *S=NULL;
16956affb6f8SHong Zhang   else {
16966affb6f8SHong Zhang     ierr = (*ts->ops->forwardgetstages)(ts,ns,S);CHKERRQ(ierr);
16976affb6f8SHong Zhang   }
16986affb6f8SHong Zhang   PetscFunctionReturn(0);
16996affb6f8SHong Zhang }
1700