xref: /petsc/src/ts/interface/sensitivity/tssen.c (revision 7adebdde874734bc3a04bfb03517869e99da1cdf)
1814e85d6SHong Zhang #include <petsc/private/tsimpl.h>        /*I "petscts.h"  I*/
2814e85d6SHong Zhang #include <petscdraw.h>
3814e85d6SHong Zhang 
4814e85d6SHong Zhang PetscLogEvent TS_AdjointStep, TS_ForwardStep;
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);
46814e85d6SHong Zhang     ierr = MatDestroy(&ts->Jacp);CHKERRQ(ierr);
47814e85d6SHong Zhang     ts->Jacp = 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;
70814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
71c9ad9de0SHong Zhang   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
72814e85d6SHong Zhang   PetscValidPointer(Amat,4);
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
81814e85d6SHong Zhang     TSSetCostIntegrand - Sets the routine for evaluating the integral term in one or more cost functions
82814e85d6SHong Zhang 
83814e85d6SHong Zhang     Logically Collective on TS
84814e85d6SHong Zhang 
85814e85d6SHong Zhang     Input Parameters:
86814e85d6SHong Zhang +   ts - the TS context obtained from TSCreate()
87814e85d6SHong Zhang .   numcost - number of gradients to be computed, this is the number of cost functions
88814e85d6SHong Zhang .   costintegral - vector that stores the integral values
89814e85d6SHong Zhang .   rf - routine for evaluating the integrand function
90c9ad9de0SHong Zhang .   drduf - function that computes the gradients of the r's with respect to u
91814e85d6SHong 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)
92814e85d6SHong Zhang .   fwd - flag indicating whether to evaluate cost integral in the forward run or the adjoint run
93814e85d6SHong Zhang -   ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
94814e85d6SHong Zhang 
95814e85d6SHong Zhang     Calling sequence of rf:
96c9ad9de0SHong Zhang $   PetscErrorCode rf(TS ts,PetscReal t,Vec U,Vec F,void *ctx);
97814e85d6SHong Zhang 
98c9ad9de0SHong Zhang     Calling sequence of drduf:
99c9ad9de0SHong Zhang $   PetscErroCode drduf(TS ts,PetscReal t,Vec U,Vec *dRdU,void *ctx);
100814e85d6SHong Zhang 
101814e85d6SHong Zhang     Calling sequence of drdpf:
102c9ad9de0SHong Zhang $   PetscErroCode drdpf(TS ts,PetscReal t,Vec U,Vec *dRdP,void *ctx);
103814e85d6SHong Zhang 
104814e85d6SHong Zhang     Level: intermediate
105814e85d6SHong Zhang 
106814e85d6SHong Zhang     Notes:
107814e85d6SHong Zhang     For optimization there is usually a single cost function (numcost = 1). For sensitivities there may be multiple cost functions
108814e85d6SHong Zhang 
109814e85d6SHong Zhang .keywords: TS, sensitivity analysis, timestep, set, quadrature, function
110814e85d6SHong Zhang 
111814e85d6SHong Zhang .seealso: TSSetRHSJacobianP(), TSGetCostGradients(), TSSetCostGradients()
112814e85d6SHong Zhang @*/
113814e85d6SHong Zhang PetscErrorCode TSSetCostIntegrand(TS ts,PetscInt numcost,Vec costintegral,PetscErrorCode (*rf)(TS,PetscReal,Vec,Vec,void*),
114c9ad9de0SHong Zhang                                                           PetscErrorCode (*drduf)(TS,PetscReal,Vec,Vec*,void*),
115814e85d6SHong Zhang                                                           PetscErrorCode (*drdpf)(TS,PetscReal,Vec,Vec*,void*),
116814e85d6SHong Zhang                                                           PetscBool fwd,void *ctx)
117814e85d6SHong Zhang {
118814e85d6SHong Zhang   PetscErrorCode ierr;
119814e85d6SHong Zhang 
120814e85d6SHong Zhang   PetscFunctionBegin;
121814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
122814e85d6SHong Zhang   if (costintegral) PetscValidHeaderSpecific(costintegral,VEC_CLASSID,3);
123814e85d6SHong 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()");
124814e85d6SHong Zhang   if (!ts->numcost) ts->numcost=numcost;
125814e85d6SHong Zhang 
126814e85d6SHong Zhang   if (costintegral) {
127814e85d6SHong Zhang     ierr = PetscObjectReference((PetscObject)costintegral);CHKERRQ(ierr);
128814e85d6SHong Zhang     ierr = VecDestroy(&ts->vec_costintegral);CHKERRQ(ierr);
129814e85d6SHong Zhang     ts->vec_costintegral = costintegral;
130814e85d6SHong Zhang   } else {
131814e85d6SHong Zhang     if (!ts->vec_costintegral) { /* Create a seq vec if user does not provide one */
132814e85d6SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,numcost,&ts->vec_costintegral);CHKERRQ(ierr);
133814e85d6SHong Zhang     } else {
134814e85d6SHong Zhang       ierr = VecSet(ts->vec_costintegral,0.0);CHKERRQ(ierr);
135814e85d6SHong Zhang     }
136814e85d6SHong Zhang   }
137814e85d6SHong Zhang   if (!ts->vec_costintegrand) {
138814e85d6SHong Zhang     ierr = VecDuplicate(ts->vec_costintegral,&ts->vec_costintegrand);CHKERRQ(ierr);
139814e85d6SHong Zhang   } else {
140814e85d6SHong Zhang     ierr = VecSet(ts->vec_costintegrand,0.0);CHKERRQ(ierr);
141814e85d6SHong Zhang   }
142814e85d6SHong Zhang   ts->costintegralfwd  = fwd; /* Evaluate the cost integral in forward run if fwd is true */
143814e85d6SHong Zhang   ts->costintegrand    = rf;
144814e85d6SHong Zhang   ts->costintegrandctx = ctx;
145c9ad9de0SHong Zhang   ts->drdufunction     = drduf;
146814e85d6SHong Zhang   ts->drdpfunction     = drdpf;
147814e85d6SHong Zhang   PetscFunctionReturn(0);
148814e85d6SHong Zhang }
149814e85d6SHong Zhang 
150b5b4017aSHong Zhang /*@C
151814e85d6SHong Zhang    TSGetCostIntegral - Returns the values of the integral term in the cost functions.
152814e85d6SHong Zhang    It is valid to call the routine after a backward run.
153814e85d6SHong Zhang 
154814e85d6SHong Zhang    Not Collective
155814e85d6SHong Zhang 
156814e85d6SHong Zhang    Input Parameter:
157814e85d6SHong Zhang .  ts - the TS context obtained from TSCreate()
158814e85d6SHong Zhang 
159814e85d6SHong Zhang    Output Parameter:
160814e85d6SHong Zhang .  v - the vector containing the integrals for each cost function
161814e85d6SHong Zhang 
162814e85d6SHong Zhang    Level: intermediate
163814e85d6SHong Zhang 
164814e85d6SHong Zhang .seealso: TSSetCostIntegrand()
165814e85d6SHong Zhang 
166814e85d6SHong Zhang .keywords: TS, sensitivity analysis
167814e85d6SHong Zhang @*/
168814e85d6SHong Zhang PetscErrorCode  TSGetCostIntegral(TS ts,Vec *v)
169814e85d6SHong Zhang {
170814e85d6SHong Zhang   PetscFunctionBegin;
171814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
172814e85d6SHong Zhang   PetscValidPointer(v,2);
173814e85d6SHong Zhang   *v = ts->vec_costintegral;
174814e85d6SHong Zhang   PetscFunctionReturn(0);
175814e85d6SHong Zhang }
176814e85d6SHong Zhang 
177b5b4017aSHong Zhang /*@C
178814e85d6SHong Zhang    TSComputeCostIntegrand - Evaluates the integral function in the cost functions.
179814e85d6SHong Zhang 
180814e85d6SHong Zhang    Input Parameters:
181814e85d6SHong Zhang +  ts - the TS context
182814e85d6SHong Zhang .  t - current time
183c9ad9de0SHong Zhang -  U - state vector, i.e. current solution
184814e85d6SHong Zhang 
185814e85d6SHong Zhang    Output Parameter:
186c9ad9de0SHong Zhang .  Q - vector of size numcost to hold the outputs
187814e85d6SHong Zhang 
188814e85d6SHong Zhang    Note:
189814e85d6SHong Zhang    Most users should not need to explicitly call this routine, as it
190814e85d6SHong Zhang    is used internally within the sensitivity analysis context.
191814e85d6SHong Zhang 
192814e85d6SHong Zhang    Level: developer
193814e85d6SHong Zhang 
194814e85d6SHong Zhang .keywords: TS, compute
195814e85d6SHong Zhang 
196814e85d6SHong Zhang .seealso: TSSetCostIntegrand()
197814e85d6SHong Zhang @*/
198c9ad9de0SHong Zhang PetscErrorCode TSComputeCostIntegrand(TS ts,PetscReal t,Vec U,Vec Q)
199814e85d6SHong Zhang {
200814e85d6SHong Zhang   PetscErrorCode ierr;
201814e85d6SHong Zhang 
202814e85d6SHong Zhang   PetscFunctionBegin;
203814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
204c9ad9de0SHong Zhang   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
205c9ad9de0SHong Zhang   PetscValidHeaderSpecific(Q,VEC_CLASSID,4);
206814e85d6SHong Zhang 
207c9ad9de0SHong Zhang   ierr = PetscLogEventBegin(TS_FunctionEval,ts,U,Q,0);CHKERRQ(ierr);
208814e85d6SHong Zhang   if (ts->costintegrand) {
209814e85d6SHong Zhang     PetscStackPush("TS user integrand in the cost function");
210c9ad9de0SHong Zhang     ierr = (*ts->costintegrand)(ts,t,U,Q,ts->costintegrandctx);CHKERRQ(ierr);
211814e85d6SHong Zhang     PetscStackPop;
212814e85d6SHong Zhang   } else {
213c9ad9de0SHong Zhang     ierr = VecZeroEntries(Q);CHKERRQ(ierr);
214814e85d6SHong Zhang   }
215814e85d6SHong Zhang 
216c9ad9de0SHong Zhang   ierr = PetscLogEventEnd(TS_FunctionEval,ts,U,Q,0);CHKERRQ(ierr);
217814e85d6SHong Zhang   PetscFunctionReturn(0);
218814e85d6SHong Zhang }
219814e85d6SHong Zhang 
220b5b4017aSHong Zhang /*@C
221c9ad9de0SHong Zhang   TSComputeDRDUFunction - Runs the user-defined DRDU function.
222814e85d6SHong Zhang 
223814e85d6SHong Zhang   Collective on TS
224814e85d6SHong Zhang 
225814e85d6SHong Zhang   Input Parameters:
226c9ad9de0SHong Zhang + ts - the TS context obtained from TSCreate()
227c9ad9de0SHong Zhang . t - current time
228c9ad9de0SHong Zhang - U - stata vector
229c9ad9de0SHong Zhang 
230c9ad9de0SHong Zhang   Output Parameters:
231b5b4017aSHong Zhang . DRDU - vector array to hold the outputs
232814e85d6SHong Zhang 
233814e85d6SHong Zhang   Notes:
234c9ad9de0SHong Zhang   TSComputeDRDUFunction() is typically used for sensitivity implementation,
235814e85d6SHong Zhang   so most users would not generally call this routine themselves.
236814e85d6SHong Zhang 
237814e85d6SHong Zhang   Level: developer
238814e85d6SHong Zhang 
239814e85d6SHong Zhang .keywords: TS, sensitivity
240814e85d6SHong Zhang .seealso: TSSetCostIntegrand()
241814e85d6SHong Zhang @*/
242c9ad9de0SHong Zhang PetscErrorCode TSComputeDRDUFunction(TS ts,PetscReal t,Vec U,Vec *DRDU)
243814e85d6SHong Zhang {
244814e85d6SHong Zhang   PetscErrorCode ierr;
245814e85d6SHong Zhang 
246814e85d6SHong Zhang   PetscFunctionBegin;
247814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
248c9ad9de0SHong Zhang   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
249814e85d6SHong Zhang 
250c9ad9de0SHong Zhang   PetscStackPush("TS user DRDU function for sensitivity analysis");
251c9ad9de0SHong Zhang   ierr = (*ts->drdufunction)(ts,t,U,DRDU,ts->costintegrandctx); CHKERRQ(ierr);
252814e85d6SHong Zhang   PetscStackPop;
253814e85d6SHong Zhang   PetscFunctionReturn(0);
254814e85d6SHong Zhang }
255814e85d6SHong Zhang 
256b5b4017aSHong Zhang /*@C
257814e85d6SHong Zhang   TSComputeDRDPFunction - Runs the user-defined DRDP function.
258814e85d6SHong Zhang 
259814e85d6SHong Zhang   Collective on TS
260814e85d6SHong Zhang 
261814e85d6SHong Zhang   Input Parameters:
262c9ad9de0SHong Zhang + ts - the TS context obtained from TSCreate()
263c9ad9de0SHong Zhang . t - current time
264c9ad9de0SHong Zhang - U - stata vector
265c9ad9de0SHong Zhang 
266c9ad9de0SHong Zhang   Output Parameters:
267b5b4017aSHong Zhang . DRDP - vector array to hold the outputs
268814e85d6SHong Zhang 
269814e85d6SHong Zhang   Notes:
270c9ad9de0SHong Zhang   TSComputeDRDPFunction() is typically used for sensitivity implementation,
271814e85d6SHong Zhang   so most users would not generally call this routine themselves.
272814e85d6SHong Zhang 
273814e85d6SHong Zhang   Level: developer
274814e85d6SHong Zhang 
275814e85d6SHong Zhang .keywords: TS, sensitivity
276814e85d6SHong Zhang .seealso: TSSetCostIntegrand()
277814e85d6SHong Zhang @*/
278c9ad9de0SHong Zhang PetscErrorCode TSComputeDRDPFunction(TS ts,PetscReal t,Vec U,Vec *DRDP)
279814e85d6SHong Zhang {
280814e85d6SHong Zhang   PetscErrorCode ierr;
281814e85d6SHong Zhang 
282814e85d6SHong Zhang   PetscFunctionBegin;
283814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
284c9ad9de0SHong Zhang   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
285814e85d6SHong Zhang 
286814e85d6SHong Zhang   PetscStackPush("TS user DRDP function for sensitivity analysis");
287c9ad9de0SHong Zhang   ierr = (*ts->drdpfunction)(ts,t,U,DRDP,ts->costintegrandctx); CHKERRQ(ierr);
288814e85d6SHong Zhang   PetscStackPop;
289814e85d6SHong Zhang   PetscFunctionReturn(0);
290814e85d6SHong Zhang }
291814e85d6SHong Zhang 
292b5b4017aSHong Zhang /*@C
293b5b4017aSHong 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.
294b5b4017aSHong Zhang 
295b5b4017aSHong Zhang   Logically Collective on TS
296b5b4017aSHong Zhang 
297b5b4017aSHong Zhang   Input Parameters:
298b5b4017aSHong Zhang + ts - TS context obtained from TSCreate()
299b5b4017aSHong Zhang . ihp1 - an array of vectors storing the result of vector-Hessian-vector product for F_UU
300b5b4017aSHong Zhang . hessianproductfunc1 - vector-Hessian-vector product function for F_UU
301b5b4017aSHong Zhang . ihp2 - an array of vectors storing the result of vector-Hessian-vector product for F_UP
302b5b4017aSHong Zhang . hessianproductfunc2 - vector-Hessian-vector product function for F_UP
303b5b4017aSHong Zhang . ihp3 - an array of vectors storing the result of vector-Hessian-vector product for F_PU
304b5b4017aSHong Zhang . hessianproductfunc3 - vector-Hessian-vector product function for F_PU
305b5b4017aSHong Zhang . ihp4 - an array of vectors storing the result of vector-Hessian-vector product for F_PP
306b5b4017aSHong Zhang . hessianproductfunc4 - vector-Hessian-vector product function for F_PP
307b5b4017aSHong Zhang 
308b5b4017aSHong Zhang   Calling sequence of ihessianproductfunc:
309b5b4017aSHong Zhang $ ihessianproductfunc (TS ts,PetscReal t,Vec U,Vec Vl,Vec Vr,Vec VHV,void *ctx);
310b5b4017aSHong Zhang +   t - current timestep
311b5b4017aSHong Zhang .   U - input vector (current ODE solution)
312b5b4017aSHong Zhang .   Vl - input vector to be left-multiplied with the Hessian
313b5b4017aSHong Zhang .   Vr - input vector to be right-multiplied with the Hessian
314b5b4017aSHong Zhang .   VHV - output vector for vector-Hessian-vector product
315b5b4017aSHong Zhang -   ctx - [optional] user-defined function context
316b5b4017aSHong Zhang 
317b5b4017aSHong Zhang   Level: intermediate
318b5b4017aSHong Zhang 
319b5b4017aSHong Zhang Note: The first Hessian function and the working array are required.
320b5b4017aSHong Zhang 
321b5b4017aSHong Zhang .keywords: TS, sensitivity
322b5b4017aSHong Zhang 
323b5b4017aSHong Zhang .seealso:
324b5b4017aSHong Zhang @*/
325b5b4017aSHong Zhang PetscErrorCode TSSetIHessianProduct(TS ts,Vec *ihp1,PetscErrorCode (*ihessianproductfunc1)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*),
326b5b4017aSHong Zhang                                           Vec *ihp2,PetscErrorCode (*ihessianproductfunc2)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*),
327b5b4017aSHong Zhang                                           Vec *ihp3,PetscErrorCode (*ihessianproductfunc3)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*),
328b5b4017aSHong Zhang                                           Vec *ihp4,PetscErrorCode (*ihessianproductfunc4)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*),
329b5b4017aSHong Zhang                                     void *ctx)
330b5b4017aSHong Zhang {
331b5b4017aSHong Zhang   PetscFunctionBegin;
332b5b4017aSHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
333b5b4017aSHong Zhang   PetscValidPointer(ihp1,2);
334b5b4017aSHong Zhang 
335b5b4017aSHong Zhang   ts->ihessianproductctx = ctx;
336b5b4017aSHong Zhang   if (ihp1) ts->vecs_fuu = ihp1;
337b5b4017aSHong Zhang   if (ihp2) ts->vecs_fup = ihp2;
338b5b4017aSHong Zhang   if (ihp3) ts->vecs_fpu = ihp3;
339b5b4017aSHong Zhang   if (ihp4) ts->vecs_fpp = ihp4;
340b5b4017aSHong Zhang   ts->ihessianproduct_fuu = ihessianproductfunc1;
341b5b4017aSHong Zhang   ts->ihessianproduct_fup = ihessianproductfunc2;
342b5b4017aSHong Zhang   ts->ihessianproduct_fpu = ihessianproductfunc3;
343b5b4017aSHong Zhang   ts->ihessianproduct_fpp = ihessianproductfunc4;
344b5b4017aSHong Zhang   PetscFunctionReturn(0);
345b5b4017aSHong Zhang }
346b5b4017aSHong Zhang 
347b5b4017aSHong Zhang /*@C
348b5b4017aSHong Zhang   TSComputeIHessianProductFunction1 - Runs the user-defined vector-Hessian-vector product function for Fuu.
349b5b4017aSHong Zhang 
350b5b4017aSHong Zhang   Collective on TS
351b5b4017aSHong Zhang 
352b5b4017aSHong Zhang   Input Parameters:
353b5b4017aSHong Zhang . ts   - The TS context obtained from TSCreate()
354b5b4017aSHong Zhang 
355b5b4017aSHong Zhang   Notes:
356b5b4017aSHong Zhang   TSComputeIHessianProductFunction1() is typically used for sensitivity implementation,
357b5b4017aSHong Zhang   so most users would not generally call this routine themselves.
358b5b4017aSHong Zhang 
359b5b4017aSHong Zhang   Level: developer
360b5b4017aSHong Zhang 
361b5b4017aSHong Zhang .keywords: TS, sensitivity
362b5b4017aSHong Zhang 
363b5b4017aSHong Zhang .seealso: TSSetIHessianProduct()
364b5b4017aSHong Zhang @*/
365b5b4017aSHong Zhang PetscErrorCode TSComputeIHessianProductFunction1(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV)
366b5b4017aSHong Zhang {
367b5b4017aSHong Zhang   PetscErrorCode ierr;
368b5b4017aSHong Zhang 
369b5b4017aSHong Zhang   PetscFunctionBegin;
370b5b4017aSHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
371b5b4017aSHong Zhang   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
372b5b4017aSHong Zhang 
373b5b4017aSHong Zhang   PetscStackPush("TS user IHessianProduct function 1 for sensitivity analysis");
374b5b4017aSHong Zhang   ierr = (*ts->ihessianproduct_fuu)(ts,t,U,Vl,Vr,VHV,ts->ihessianproductctx);CHKERRQ(ierr);
375b5b4017aSHong Zhang   PetscStackPop;
376b5b4017aSHong Zhang   PetscFunctionReturn(0);
377b5b4017aSHong Zhang }
378b5b4017aSHong Zhang 
379b5b4017aSHong Zhang /*@C
380b5b4017aSHong Zhang   TSComputeIHessianProductFunction2 - Runs the user-defined vector-Hessian-vector product function for Fup.
381b5b4017aSHong Zhang 
382b5b4017aSHong Zhang   Collective on TS
383b5b4017aSHong Zhang 
384b5b4017aSHong Zhang   Input Parameters:
385b5b4017aSHong Zhang . ts   - The TS context obtained from TSCreate()
386b5b4017aSHong Zhang 
387b5b4017aSHong Zhang   Notes:
388b5b4017aSHong Zhang   TSComputeIHessianProductFunction2() is typically used for sensitivity implementation,
389b5b4017aSHong Zhang   so most users would not generally call this routine themselves.
390b5b4017aSHong Zhang 
391b5b4017aSHong Zhang   Level: developer
392b5b4017aSHong Zhang 
393b5b4017aSHong Zhang .keywords: TS, sensitivity
394b5b4017aSHong Zhang 
395b5b4017aSHong Zhang .seealso: TSSetIHessianProduct()
396b5b4017aSHong Zhang @*/
397b5b4017aSHong Zhang PetscErrorCode TSComputeIHessianProductFunction2(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV)
398b5b4017aSHong Zhang {
399b5b4017aSHong Zhang   PetscErrorCode ierr;
400b5b4017aSHong Zhang 
401b5b4017aSHong Zhang   PetscFunctionBegin;
402b5b4017aSHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
403b5b4017aSHong Zhang   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
404b5b4017aSHong Zhang 
405b5b4017aSHong Zhang   PetscStackPush("TS user IHessianProduct function 2 for sensitivity analysis");
406b5b4017aSHong Zhang   ierr = (*ts->ihessianproduct_fup)(ts,t,U,Vl,Vr,VHV,ts->ihessianproductctx);CHKERRQ(ierr);
407b5b4017aSHong Zhang   PetscStackPop;
408b5b4017aSHong Zhang   PetscFunctionReturn(0);
409b5b4017aSHong Zhang }
410b5b4017aSHong Zhang 
411b5b4017aSHong Zhang /*@C
412b5b4017aSHong Zhang   TSComputeIHessianProductFunction3 - Runs the user-defined vector-Hessian-vector product function for Fpu.
413b5b4017aSHong Zhang 
414b5b4017aSHong Zhang   Collective on TS
415b5b4017aSHong Zhang 
416b5b4017aSHong Zhang   Input Parameters:
417b5b4017aSHong Zhang . ts   - The TS context obtained from TSCreate()
418b5b4017aSHong Zhang 
419b5b4017aSHong Zhang   Notes:
420b5b4017aSHong Zhang   TSComputeIHessianProductFunction3() is typically used for sensitivity implementation,
421b5b4017aSHong Zhang   so most users would not generally call this routine themselves.
422b5b4017aSHong Zhang 
423b5b4017aSHong Zhang   Level: developer
424b5b4017aSHong Zhang 
425b5b4017aSHong Zhang .keywords: TS, sensitivity
426b5b4017aSHong Zhang 
427b5b4017aSHong Zhang .seealso: TSSetIHessianProduct()
428b5b4017aSHong Zhang @*/
429b5b4017aSHong Zhang PetscErrorCode TSComputeIHessianProductFunction3(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV)
430b5b4017aSHong Zhang {
431b5b4017aSHong Zhang   PetscErrorCode ierr;
432b5b4017aSHong Zhang 
433b5b4017aSHong Zhang   PetscFunctionBegin;
434b5b4017aSHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
435b5b4017aSHong Zhang   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
436b5b4017aSHong Zhang 
437b5b4017aSHong Zhang   PetscStackPush("TS user IHessianProduct function 3 for sensitivity analysis");
438b5b4017aSHong Zhang   ierr = (*ts->ihessianproduct_fpu)(ts,t,U,Vl,Vr,VHV,ts->ihessianproductctx);CHKERRQ(ierr);
439b5b4017aSHong Zhang   PetscStackPop;
440b5b4017aSHong Zhang   PetscFunctionReturn(0);
441b5b4017aSHong Zhang }
442b5b4017aSHong Zhang 
443b5b4017aSHong Zhang /*@C
444b5b4017aSHong Zhang   TSComputeIHessianProductFunction4 - Runs the user-defined vector-Hessian-vector product function for Fpp.
445b5b4017aSHong Zhang 
446b5b4017aSHong Zhang   Collective on TS
447b5b4017aSHong Zhang 
448b5b4017aSHong Zhang   Input Parameters:
449b5b4017aSHong Zhang . ts   - The TS context obtained from TSCreate()
450b5b4017aSHong Zhang 
451b5b4017aSHong Zhang   Notes:
452b5b4017aSHong Zhang   TSComputeIHessianProductFunction4() is typically used for sensitivity implementation,
453b5b4017aSHong Zhang   so most users would not generally call this routine themselves.
454b5b4017aSHong Zhang 
455b5b4017aSHong Zhang   Level: developer
456b5b4017aSHong Zhang 
457b5b4017aSHong Zhang .keywords: TS, sensitivity
458b5b4017aSHong Zhang 
459b5b4017aSHong Zhang .seealso: TSSetIHessianProduct()
460b5b4017aSHong Zhang @*/
461b5b4017aSHong Zhang PetscErrorCode TSComputeIHessianProductFunction4(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV)
462b5b4017aSHong Zhang {
463b5b4017aSHong Zhang   PetscErrorCode ierr;
464b5b4017aSHong Zhang 
465b5b4017aSHong Zhang   PetscFunctionBegin;
466b5b4017aSHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
467b5b4017aSHong Zhang   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
468b5b4017aSHong Zhang 
469b5b4017aSHong Zhang   PetscStackPush("TS user IHessianProduct function 3 for sensitivity analysis");
470b5b4017aSHong Zhang   ierr = (*ts->ihessianproduct_fpp)(ts,t,U,Vl,Vr,VHV,ts->ihessianproductctx);CHKERRQ(ierr);
471b5b4017aSHong Zhang   PetscStackPop;
472b5b4017aSHong Zhang   PetscFunctionReturn(0);
473b5b4017aSHong Zhang }
474b5b4017aSHong Zhang 
475814e85d6SHong Zhang /* --------------------------- Adjoint sensitivity ---------------------------*/
476814e85d6SHong Zhang 
477814e85d6SHong Zhang /*@
478814e85d6SHong 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
479814e85d6SHong Zhang       for use by the TSAdjoint routines.
480814e85d6SHong Zhang 
481814e85d6SHong Zhang    Logically Collective on TS and Vec
482814e85d6SHong Zhang 
483814e85d6SHong Zhang    Input Parameters:
484814e85d6SHong Zhang +  ts - the TS context obtained from TSCreate()
485814e85d6SHong 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
486814e85d6SHong Zhang -  mu - gradients with respect to the parameters, the number of entries in these vectors is the same as the number of parameters
487814e85d6SHong Zhang 
488814e85d6SHong Zhang    Level: beginner
489814e85d6SHong Zhang 
490814e85d6SHong Zhang    Notes:
491814e85d6SHong Zhang     the entries in these vectors must be correctly initialized with the values lamda_i = df/dy|finaltime  mu_i = df/dp|finaltime
492814e85d6SHong Zhang 
493814e85d6SHong Zhang    After TSAdjointSolve() is called the lamba and the mu contain the computed sensitivities
494814e85d6SHong Zhang 
495814e85d6SHong Zhang .keywords: TS, timestep, set, sensitivity, initial values
496b5b4017aSHong Zhang 
497b5b4017aSHong Zhang .seealso TSGetCostGradients()
498814e85d6SHong Zhang @*/
499814e85d6SHong Zhang PetscErrorCode TSSetCostGradients(TS ts,PetscInt numcost,Vec *lambda,Vec *mu)
500814e85d6SHong Zhang {
501814e85d6SHong Zhang   PetscFunctionBegin;
502814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
503814e85d6SHong Zhang   PetscValidPointer(lambda,2);
504814e85d6SHong Zhang   ts->vecs_sensi  = lambda;
505814e85d6SHong Zhang   ts->vecs_sensip = mu;
506814e85d6SHong 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");
507814e85d6SHong Zhang   ts->numcost  = numcost;
508814e85d6SHong Zhang   PetscFunctionReturn(0);
509814e85d6SHong Zhang }
510814e85d6SHong Zhang 
511814e85d6SHong Zhang /*@
512814e85d6SHong Zhang    TSGetCostGradients - Returns the gradients from the TSAdjointSolve()
513814e85d6SHong Zhang 
514814e85d6SHong Zhang    Not Collective, but Vec returned is parallel if TS is parallel
515814e85d6SHong Zhang 
516814e85d6SHong Zhang    Input Parameter:
517814e85d6SHong Zhang .  ts - the TS context obtained from TSCreate()
518814e85d6SHong Zhang 
519814e85d6SHong Zhang    Output Parameter:
520814e85d6SHong Zhang +  lambda - vectors containing the gradients of the cost functions with respect to the ODE/DAE solution variables
521814e85d6SHong Zhang -  mu - vectors containing the gradients of the cost functions with respect to the problem parameters
522814e85d6SHong Zhang 
523814e85d6SHong Zhang    Level: intermediate
524814e85d6SHong Zhang 
525814e85d6SHong Zhang .keywords: TS, timestep, get, sensitivity
526b5b4017aSHong Zhang 
527b5b4017aSHong Zhang .seealso: TSSetCostGradients()
528814e85d6SHong Zhang @*/
529814e85d6SHong Zhang PetscErrorCode TSGetCostGradients(TS ts,PetscInt *numcost,Vec **lambda,Vec **mu)
530814e85d6SHong Zhang {
531814e85d6SHong Zhang   PetscFunctionBegin;
532814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
533814e85d6SHong Zhang   if (numcost) *numcost = ts->numcost;
534814e85d6SHong Zhang   if (lambda)  *lambda  = ts->vecs_sensi;
535814e85d6SHong Zhang   if (mu)      *mu      = ts->vecs_sensip;
536814e85d6SHong Zhang   PetscFunctionReturn(0);
537814e85d6SHong Zhang }
538814e85d6SHong Zhang 
539814e85d6SHong Zhang /*@
540b5b4017aSHong 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
541b5b4017aSHong Zhang       for use by the TSAdjoint routines.
542b5b4017aSHong Zhang 
543b5b4017aSHong Zhang    Logically Collective on TS and Vec
544b5b4017aSHong Zhang 
545b5b4017aSHong Zhang    Input Parameters:
546b5b4017aSHong Zhang +  ts - the TS context obtained from TSCreate()
547b5b4017aSHong Zhang .  numcost - number of cost functions
548b5b4017aSHong 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
549b5b4017aSHong 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
550b5b4017aSHong Zhang -  dir - the direction vector that are multiplied with the Hessian of the cost functions
551b5b4017aSHong Zhang 
552b5b4017aSHong Zhang    Level: beginner
553b5b4017aSHong Zhang 
554b5b4017aSHong Zhang    Notes: Hessian of the cost function is completely different from Hessian of the ODE/DAE system
555b5b4017aSHong Zhang 
556b5b4017aSHong Zhang    For second-order adjoint, one needs to call this function and then TSAdjointInitializeForward() before TSSolve().
557b5b4017aSHong Zhang 
558b5b4017aSHong 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.
559b5b4017aSHong Zhang 
560b5b4017aSHong Zhang .keywords: TS, sensitivity, second-order adjoint
561b5b4017aSHong Zhang 
562b5b4017aSHong Zhang .seealso: TSAdjointInitializeForward()
563b5b4017aSHong Zhang @*/
564b5b4017aSHong Zhang PetscErrorCode TSSetCostHessianProducts(TS ts,PetscInt numcost,Vec *lambda2,Vec *mu2,Vec dir)
565b5b4017aSHong Zhang {
566b5b4017aSHong Zhang   PetscFunctionBegin;
567b5b4017aSHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
568b5b4017aSHong Zhang   PetscValidPointer(lambda2,2);
569b5b4017aSHong 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");
570b5b4017aSHong Zhang   ts->numcost       = numcost;
571b5b4017aSHong Zhang   ts->vecs_sensi2   = lambda2;
572b5b4017aSHong Zhang   ts->vecs_sensip2  = mu2;
573b5b4017aSHong Zhang   ts->vec_dir       = dir;
574b5b4017aSHong Zhang   PetscFunctionReturn(0);
575b5b4017aSHong Zhang }
576b5b4017aSHong Zhang 
577b5b4017aSHong Zhang /*@
578b5b4017aSHong Zhang    TSGetCostHessianProducts - Returns the gradients from the TSAdjointSolve()
579b5b4017aSHong Zhang 
580b5b4017aSHong Zhang    Not Collective, but Vec returned is parallel if TS is parallel
581b5b4017aSHong Zhang 
582b5b4017aSHong Zhang    Input Parameter:
583b5b4017aSHong Zhang .  ts - the TS context obtained from TSCreate()
584b5b4017aSHong Zhang 
585b5b4017aSHong Zhang    Output Parameter:
586b5b4017aSHong Zhang +  numcost - number of cost functions
587b5b4017aSHong 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
588b5b4017aSHong 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
589b5b4017aSHong Zhang -  dir - the direction vector that are multiplied with the Hessian of the cost functions
590b5b4017aSHong Zhang 
591b5b4017aSHong Zhang    Level: intermediate
592b5b4017aSHong Zhang 
593b5b4017aSHong Zhang .keywords: TS, sensitivity, second-order adjoint
594b5b4017aSHong Zhang 
595b5b4017aSHong Zhang .seealso: TSSetCostHessianProducts()
596b5b4017aSHong Zhang @*/
597b5b4017aSHong Zhang PetscErrorCode TSGetCostHessianProducts(TS ts,PetscInt *numcost,Vec **lambda2,Vec **mu2, Vec *dir)
598b5b4017aSHong Zhang {
599b5b4017aSHong Zhang   PetscFunctionBegin;
600b5b4017aSHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
601b5b4017aSHong Zhang   if (numcost) *numcost = ts->numcost;
602b5b4017aSHong Zhang   if (lambda2) *lambda2 = ts->vecs_sensi2;
603b5b4017aSHong Zhang   if (mu2)     *mu2     = ts->vecs_sensip2;
604b5b4017aSHong Zhang   if (dir)     *dir     = ts->vec_dir;
605b5b4017aSHong Zhang   PetscFunctionReturn(0);
606b5b4017aSHong Zhang }
607b5b4017aSHong Zhang 
608b5b4017aSHong Zhang /*@
609b5b4017aSHong Zhang   TSAdjointInitializeForward - Trigger the tangent linear solver and initialize the forward sensitivities
610b5b4017aSHong Zhang 
611b5b4017aSHong Zhang   Logically Collective on TS and Mat
612b5b4017aSHong Zhang 
613b5b4017aSHong Zhang   Input Parameters:
614b5b4017aSHong Zhang +  ts - the TS context obtained from TSCreate()
615b5b4017aSHong Zhang -  didp - the derivative of initial values w.r.t. parameters
616b5b4017aSHong Zhang 
617b5b4017aSHong Zhang   Level: intermediate
618b5b4017aSHong Zhang 
619b5b4017aSHong 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.
620b5b4017aSHong Zhang 
621b5b4017aSHong Zhang .keywords: TS, sensitivity, second-order adjoint
622b5b4017aSHong Zhang 
623b5b4017aSHong Zhang .seealso: TSSetCostHessianProducts()
624b5b4017aSHong Zhang @*/
625b5b4017aSHong Zhang PetscErrorCode TSAdjointInitializeForward(TS ts,Mat didp)
626b5b4017aSHong Zhang {
627b5b4017aSHong Zhang   PetscErrorCode ierr;
628b5b4017aSHong Zhang 
629b5b4017aSHong Zhang   PetscFunctionBegin;
630b5b4017aSHong Zhang   ts->forward_solve = PETSC_TRUE; /* turn on tangent linear mode */
631b5b4017aSHong Zhang   if (!ts->vecs_sensi2) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must call TSSetCostHessianProducts() first");
632b5b4017aSHong Zhang   if (ts->vecs_sensip2 && !didp) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"The fourth argument is not NULL, indicating parametric sensitivities are desired, so the dIdP matrix must be provided"); /* check conflicted settings */
633b5b4017aSHong Zhang   ierr = TSForwardSetInitialSensitivities(ts,didp);CHKERRQ(ierr); /* if didp is NULL, identity matrix is assumed */
634b5b4017aSHong Zhang   PetscFunctionReturn(0);
635b5b4017aSHong Zhang }
636b5b4017aSHong Zhang 
637b5b4017aSHong Zhang /*@
638814e85d6SHong Zhang    TSAdjointSetUp - Sets up the internal data structures for the later use
639814e85d6SHong Zhang    of an adjoint solver
640814e85d6SHong Zhang 
641814e85d6SHong Zhang    Collective on TS
642814e85d6SHong Zhang 
643814e85d6SHong Zhang    Input Parameter:
644814e85d6SHong Zhang .  ts - the TS context obtained from TSCreate()
645814e85d6SHong Zhang 
646814e85d6SHong Zhang    Level: advanced
647814e85d6SHong Zhang 
648814e85d6SHong Zhang .keywords: TS, timestep, setup
649814e85d6SHong Zhang 
650814e85d6SHong Zhang .seealso: TSCreate(), TSAdjointStep(), TSSetCostGradients()
651814e85d6SHong Zhang @*/
652814e85d6SHong Zhang PetscErrorCode TSAdjointSetUp(TS ts)
653814e85d6SHong Zhang {
654814e85d6SHong Zhang   PetscErrorCode ierr;
655814e85d6SHong Zhang 
656814e85d6SHong Zhang   PetscFunctionBegin;
657814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
658814e85d6SHong Zhang   if (ts->adjointsetupcalled) PetscFunctionReturn(0);
659814e85d6SHong Zhang   if (!ts->vecs_sensi) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetCostGradients() first");
660814e85d6SHong Zhang   if (ts->vecs_sensip && !ts->Jacp) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONGSTATE,"Must call TSAdjointSetRHSJacobian() first");
661814e85d6SHong Zhang 
662814e85d6SHong Zhang   if (ts->vec_costintegral) { /* if there is integral in the cost function */
663c9ad9de0SHong Zhang     ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&ts->vecs_drdu);CHKERRQ(ierr);
664814e85d6SHong Zhang     if (ts->vecs_sensip){
665814e85d6SHong Zhang       ierr = VecDuplicateVecs(ts->vecs_sensip[0],ts->numcost,&ts->vecs_drdp);CHKERRQ(ierr);
666814e85d6SHong Zhang     }
667814e85d6SHong Zhang   }
668814e85d6SHong Zhang 
669814e85d6SHong Zhang   if (ts->ops->adjointsetup) {
670814e85d6SHong Zhang     ierr = (*ts->ops->adjointsetup)(ts);CHKERRQ(ierr);
671814e85d6SHong Zhang   }
672814e85d6SHong Zhang   ts->adjointsetupcalled = PETSC_TRUE;
673814e85d6SHong Zhang   PetscFunctionReturn(0);
674814e85d6SHong Zhang }
675814e85d6SHong Zhang 
676814e85d6SHong Zhang /*@
677ece46509SHong Zhang    TSAdjointReset - Resets a TSAdjoint context and removes any allocated Vecs and Mats.
678ece46509SHong Zhang 
679ece46509SHong Zhang    Collective on TS
680ece46509SHong Zhang 
681ece46509SHong Zhang    Input Parameter:
682ece46509SHong Zhang .  ts - the TS context obtained from TSCreate()
683ece46509SHong Zhang 
684ece46509SHong Zhang    Level: beginner
685ece46509SHong Zhang 
686ece46509SHong Zhang .keywords: TS, timestep, reset
687ece46509SHong Zhang 
688ece46509SHong Zhang .seealso: TSCreate(), TSAdjointSetup(), TSADestroy()
689ece46509SHong Zhang @*/
690ece46509SHong Zhang PetscErrorCode TSAdjointReset(TS ts)
691ece46509SHong Zhang {
692ece46509SHong Zhang   PetscErrorCode ierr;
693ece46509SHong Zhang 
694ece46509SHong Zhang   PetscFunctionBegin;
695ece46509SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
696ece46509SHong Zhang   if (ts->ops->adjointreset) {
697ece46509SHong Zhang     ierr = (*ts->ops->adjointreset)(ts);CHKERRQ(ierr);
698ece46509SHong Zhang   }
699ece46509SHong Zhang   ts->vecs_sensi         = NULL;
700ece46509SHong Zhang   ts->vecs_sensip        = NULL;
701ece46509SHong Zhang   ts->vecs_sensi2        = NULL;
702ece46509SHong Zhang   ts->vecs_sensip2       = NULL;
703ece46509SHong Zhang   ts->vec_dir            = NULL;
704ece46509SHong Zhang   ts->adjointsetupcalled = PETSC_FALSE;
705ece46509SHong Zhang   PetscFunctionReturn(0);
706ece46509SHong Zhang }
707ece46509SHong Zhang 
708ece46509SHong Zhang /*@
709814e85d6SHong Zhang    TSAdjointSetSteps - Sets the number of steps the adjoint solver should take backward in time
710814e85d6SHong Zhang 
711814e85d6SHong Zhang    Logically Collective on TS
712814e85d6SHong Zhang 
713814e85d6SHong Zhang    Input Parameters:
714814e85d6SHong Zhang +  ts - the TS context obtained from TSCreate()
715814e85d6SHong Zhang .  steps - number of steps to use
716814e85d6SHong Zhang 
717814e85d6SHong Zhang    Level: intermediate
718814e85d6SHong Zhang 
719814e85d6SHong Zhang    Notes:
720814e85d6SHong Zhang     Normally one does not call this and TSAdjointSolve() integrates back to the original timestep. One can call this
721814e85d6SHong Zhang           so as to integrate back to less than the original timestep
722814e85d6SHong Zhang 
723814e85d6SHong Zhang .keywords: TS, timestep, set, maximum, iterations
724814e85d6SHong Zhang 
725814e85d6SHong Zhang .seealso: TSSetExactFinalTime()
726814e85d6SHong Zhang @*/
727814e85d6SHong Zhang PetscErrorCode TSAdjointSetSteps(TS ts,PetscInt steps)
728814e85d6SHong Zhang {
729814e85d6SHong Zhang   PetscFunctionBegin;
730814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
731814e85d6SHong Zhang   PetscValidLogicalCollectiveInt(ts,steps,2);
732814e85d6SHong Zhang   if (steps < 0) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Cannot step back a negative number of steps");
733814e85d6SHong Zhang   if (steps > ts->steps) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Cannot step back more than the total number of forward steps");
734814e85d6SHong Zhang   ts->adjoint_max_steps = steps;
735814e85d6SHong Zhang   PetscFunctionReturn(0);
736814e85d6SHong Zhang }
737814e85d6SHong Zhang 
738814e85d6SHong Zhang /*@C
739814e85d6SHong Zhang   TSAdjointSetRHSJacobian - Deprecated, use TSSetRHSJacobianP()
740814e85d6SHong Zhang 
741814e85d6SHong Zhang   Level: deprecated
742814e85d6SHong Zhang 
743814e85d6SHong Zhang @*/
744814e85d6SHong Zhang PetscErrorCode TSAdjointSetRHSJacobian(TS ts,Mat Amat,PetscErrorCode (*func)(TS,PetscReal,Vec,Mat,void*),void *ctx)
745814e85d6SHong Zhang {
746814e85d6SHong Zhang   PetscErrorCode ierr;
747814e85d6SHong Zhang 
748814e85d6SHong Zhang   PetscFunctionBegin;
749814e85d6SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
750814e85d6SHong Zhang   PetscValidHeaderSpecific(Amat,MAT_CLASSID,2);
751814e85d6SHong Zhang 
752814e85d6SHong Zhang   ts->rhsjacobianp    = func;
753814e85d6SHong Zhang   ts->rhsjacobianpctx = ctx;
754814e85d6SHong Zhang   if(Amat) {
755814e85d6SHong Zhang     ierr = PetscObjectReference((PetscObject)Amat);CHKERRQ(ierr);
756814e85d6SHong Zhang     ierr = MatDestroy(&ts->Jacp);CHKERRQ(ierr);
757814e85d6SHong Zhang     ts->Jacp = Amat;
758814e85d6SHong Zhang   }
759814e85d6SHong Zhang   PetscFunctionReturn(0);
760814e85d6SHong Zhang }
761814e85d6SHong Zhang 
762814e85d6SHong Zhang /*@C
763814e85d6SHong Zhang   TSAdjointComputeRHSJacobian - Deprecated, use TSComputeRHSJacobianP()
764814e85d6SHong Zhang 
765814e85d6SHong Zhang   Level: deprecated
766814e85d6SHong Zhang 
767814e85d6SHong Zhang @*/
768c9ad9de0SHong Zhang PetscErrorCode TSAdjointComputeRHSJacobian(TS ts,PetscReal t,Vec U,Mat Amat)
769814e85d6SHong Zhang {
770814e85d6SHong Zhang   PetscErrorCode ierr;
771814e85d6SHong Zhang 
772814e85d6SHong Zhang   PetscFunctionBegin;
773814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
774c9ad9de0SHong Zhang   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
775814e85d6SHong Zhang   PetscValidPointer(Amat,4);
776814e85d6SHong Zhang 
777814e85d6SHong Zhang   PetscStackPush("TS user JacobianP function for sensitivity analysis");
778c9ad9de0SHong Zhang   ierr = (*ts->rhsjacobianp)(ts,t,U,Amat,ts->rhsjacobianpctx); CHKERRQ(ierr);
779814e85d6SHong Zhang   PetscStackPop;
780814e85d6SHong Zhang   PetscFunctionReturn(0);
781814e85d6SHong Zhang }
782814e85d6SHong Zhang 
783814e85d6SHong Zhang /*@
784c9ad9de0SHong Zhang   TSAdjointComputeDRDYFunction - Deprecated, use TSComputeDRDUFunction()
785814e85d6SHong Zhang 
786814e85d6SHong Zhang   Level: deprecated
787814e85d6SHong Zhang 
788814e85d6SHong Zhang @*/
789c9ad9de0SHong Zhang PetscErrorCode TSAdjointComputeDRDYFunction(TS ts,PetscReal t,Vec U,Vec *DRDU)
790814e85d6SHong Zhang {
791814e85d6SHong Zhang   PetscErrorCode ierr;
792814e85d6SHong Zhang 
793814e85d6SHong Zhang   PetscFunctionBegin;
794814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
795c9ad9de0SHong Zhang   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
796814e85d6SHong Zhang 
797814e85d6SHong Zhang   PetscStackPush("TS user DRDY function for sensitivity analysis");
798c9ad9de0SHong Zhang   ierr = (*ts->drdufunction)(ts,t,U,DRDU,ts->costintegrandctx); CHKERRQ(ierr);
799814e85d6SHong Zhang   PetscStackPop;
800814e85d6SHong Zhang   PetscFunctionReturn(0);
801814e85d6SHong Zhang }
802814e85d6SHong Zhang 
803814e85d6SHong Zhang /*@
804814e85d6SHong Zhang   TSAdjointComputeDRDPFunction - Deprecated, use TSComputeDRDPFunction()
805814e85d6SHong Zhang 
806814e85d6SHong Zhang   Level: deprecated
807814e85d6SHong Zhang 
808814e85d6SHong Zhang @*/
809c9ad9de0SHong Zhang PetscErrorCode TSAdjointComputeDRDPFunction(TS ts,PetscReal t,Vec U,Vec *DRDP)
810814e85d6SHong Zhang {
811814e85d6SHong Zhang   PetscErrorCode ierr;
812814e85d6SHong Zhang 
813814e85d6SHong Zhang   PetscFunctionBegin;
814814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
815c9ad9de0SHong Zhang   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
816814e85d6SHong Zhang 
817814e85d6SHong Zhang   PetscStackPush("TS user DRDP function for sensitivity analysis");
818c9ad9de0SHong Zhang   ierr = (*ts->drdpfunction)(ts,t,U,DRDP,ts->costintegrandctx); CHKERRQ(ierr);
819814e85d6SHong Zhang   PetscStackPop;
820814e85d6SHong Zhang   PetscFunctionReturn(0);
821814e85d6SHong Zhang }
822814e85d6SHong Zhang 
823814e85d6SHong Zhang /*@C
824814e85d6SHong Zhang    TSAdjointMonitorSensi - monitors the first lambda sensitivity
825814e85d6SHong Zhang 
826814e85d6SHong Zhang    Level: intermediate
827814e85d6SHong Zhang 
828814e85d6SHong Zhang .keywords: TS, set, monitor
829814e85d6SHong Zhang 
830814e85d6SHong Zhang .seealso: TSAdjointMonitorSet()
831814e85d6SHong Zhang @*/
832814e85d6SHong Zhang PetscErrorCode TSAdjointMonitorSensi(TS ts,PetscInt step,PetscReal ptime,Vec v,PetscInt numcost,Vec *lambda,Vec *mu,PetscViewerAndFormat *vf)
833814e85d6SHong Zhang {
834814e85d6SHong Zhang   PetscErrorCode ierr;
835814e85d6SHong Zhang   PetscViewer    viewer = vf->viewer;
836814e85d6SHong Zhang 
837814e85d6SHong Zhang   PetscFunctionBegin;
838814e85d6SHong Zhang   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4);
839814e85d6SHong Zhang   ierr = PetscViewerPushFormat(viewer,vf->format);CHKERRQ(ierr);
840814e85d6SHong Zhang   ierr = VecView(lambda[0],viewer);CHKERRQ(ierr);
841814e85d6SHong Zhang   ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
842814e85d6SHong Zhang   PetscFunctionReturn(0);
843814e85d6SHong Zhang }
844814e85d6SHong Zhang 
845814e85d6SHong Zhang /*@C
846814e85d6SHong Zhang    TSAdjointMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user
847814e85d6SHong Zhang 
848814e85d6SHong Zhang    Collective on TS
849814e85d6SHong Zhang 
850814e85d6SHong Zhang    Input Parameters:
851814e85d6SHong Zhang +  ts - TS object you wish to monitor
852814e85d6SHong Zhang .  name - the monitor type one is seeking
853814e85d6SHong Zhang .  help - message indicating what monitoring is done
854814e85d6SHong Zhang .  manual - manual page for the monitor
855814e85d6SHong Zhang .  monitor - the monitor function
856814e85d6SHong 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
857814e85d6SHong Zhang 
858814e85d6SHong Zhang    Level: developer
859814e85d6SHong Zhang 
860814e85d6SHong Zhang .seealso: PetscOptionsGetViewer(), PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(),
861814e85d6SHong Zhang           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
862814e85d6SHong Zhang           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
863814e85d6SHong Zhang           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
864814e85d6SHong Zhang           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
865814e85d6SHong Zhang           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
866814e85d6SHong Zhang           PetscOptionsFList(), PetscOptionsEList()
867814e85d6SHong Zhang @*/
868814e85d6SHong 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*))
869814e85d6SHong Zhang {
870814e85d6SHong Zhang   PetscErrorCode    ierr;
871814e85d6SHong Zhang   PetscViewer       viewer;
872814e85d6SHong Zhang   PetscViewerFormat format;
873814e85d6SHong Zhang   PetscBool         flg;
874814e85d6SHong Zhang 
875814e85d6SHong Zhang   PetscFunctionBegin;
87616413a6aSBarry Smith   ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)ts),((PetscObject) ts)->options,((PetscObject)ts)->prefix,name,&viewer,&format,&flg);CHKERRQ(ierr);
877814e85d6SHong Zhang   if (flg) {
878814e85d6SHong Zhang     PetscViewerAndFormat *vf;
879814e85d6SHong Zhang     ierr = PetscViewerAndFormatCreate(viewer,format,&vf);CHKERRQ(ierr);
880814e85d6SHong Zhang     ierr = PetscObjectDereference((PetscObject)viewer);CHKERRQ(ierr);
881814e85d6SHong Zhang     if (monitorsetup) {
882814e85d6SHong Zhang       ierr = (*monitorsetup)(ts,vf);CHKERRQ(ierr);
883814e85d6SHong Zhang     }
884814e85d6SHong Zhang     ierr = TSAdjointMonitorSet(ts,(PetscErrorCode (*)(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*))monitor,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy);CHKERRQ(ierr);
885814e85d6SHong Zhang   }
886814e85d6SHong Zhang   PetscFunctionReturn(0);
887814e85d6SHong Zhang }
888814e85d6SHong Zhang 
889814e85d6SHong Zhang /*@C
890814e85d6SHong Zhang    TSAdjointMonitorSet - Sets an ADDITIONAL function that is to be used at every
891814e85d6SHong Zhang    timestep to display the iteration's  progress.
892814e85d6SHong Zhang 
893814e85d6SHong Zhang    Logically Collective on TS
894814e85d6SHong Zhang 
895814e85d6SHong Zhang    Input Parameters:
896814e85d6SHong Zhang +  ts - the TS context obtained from TSCreate()
897814e85d6SHong Zhang .  adjointmonitor - monitoring routine
898814e85d6SHong Zhang .  adjointmctx - [optional] user-defined context for private data for the
899814e85d6SHong Zhang              monitor routine (use NULL if no context is desired)
900814e85d6SHong Zhang -  adjointmonitordestroy - [optional] routine that frees monitor context
901814e85d6SHong Zhang           (may be NULL)
902814e85d6SHong Zhang 
903814e85d6SHong Zhang    Calling sequence of monitor:
904814e85d6SHong Zhang $    int adjointmonitor(TS ts,PetscInt steps,PetscReal time,Vec u,PetscInt numcost,Vec *lambda, Vec *mu,void *adjointmctx)
905814e85d6SHong Zhang 
906814e85d6SHong Zhang +    ts - the TS context
907814e85d6SHong 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
908814e85d6SHong Zhang                                been interpolated to)
909814e85d6SHong Zhang .    time - current time
910814e85d6SHong Zhang .    u - current iterate
911814e85d6SHong Zhang .    numcost - number of cost functionos
912814e85d6SHong Zhang .    lambda - sensitivities to initial conditions
913814e85d6SHong Zhang .    mu - sensitivities to parameters
914814e85d6SHong Zhang -    adjointmctx - [optional] adjoint monitoring context
915814e85d6SHong Zhang 
916814e85d6SHong Zhang    Notes:
917814e85d6SHong Zhang    This routine adds an additional monitor to the list of monitors that
918814e85d6SHong Zhang    already has been loaded.
919814e85d6SHong Zhang 
920814e85d6SHong Zhang    Fortran Notes:
921814e85d6SHong Zhang     Only a single monitor function can be set for each TS object
922814e85d6SHong Zhang 
923814e85d6SHong Zhang    Level: intermediate
924814e85d6SHong Zhang 
925814e85d6SHong Zhang .keywords: TS, timestep, set, adjoint, monitor
926814e85d6SHong Zhang 
927814e85d6SHong Zhang .seealso: TSAdjointMonitorCancel()
928814e85d6SHong Zhang @*/
929814e85d6SHong Zhang PetscErrorCode TSAdjointMonitorSet(TS ts,PetscErrorCode (*adjointmonitor)(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*),void *adjointmctx,PetscErrorCode (*adjointmdestroy)(void**))
930814e85d6SHong Zhang {
931814e85d6SHong Zhang   PetscErrorCode ierr;
932814e85d6SHong Zhang   PetscInt       i;
933814e85d6SHong Zhang   PetscBool      identical;
934814e85d6SHong Zhang 
935814e85d6SHong Zhang   PetscFunctionBegin;
936814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
937814e85d6SHong Zhang   for (i=0; i<ts->numbermonitors;i++) {
938814e85d6SHong Zhang     ierr = PetscMonitorCompare((PetscErrorCode (*)(void))adjointmonitor,adjointmctx,adjointmdestroy,(PetscErrorCode (*)(void))ts->adjointmonitor[i],ts->adjointmonitorcontext[i],ts->adjointmonitordestroy[i],&identical);CHKERRQ(ierr);
939814e85d6SHong Zhang     if (identical) PetscFunctionReturn(0);
940814e85d6SHong Zhang   }
941814e85d6SHong Zhang   if (ts->numberadjointmonitors >= MAXTSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many adjoint monitors set");
942814e85d6SHong Zhang   ts->adjointmonitor[ts->numberadjointmonitors]          = adjointmonitor;
943814e85d6SHong Zhang   ts->adjointmonitordestroy[ts->numberadjointmonitors]   = adjointmdestroy;
944814e85d6SHong Zhang   ts->adjointmonitorcontext[ts->numberadjointmonitors++] = (void*)adjointmctx;
945814e85d6SHong Zhang   PetscFunctionReturn(0);
946814e85d6SHong Zhang }
947814e85d6SHong Zhang 
948814e85d6SHong Zhang /*@C
949814e85d6SHong Zhang    TSAdjointMonitorCancel - Clears all the adjoint monitors that have been set on a time-step object.
950814e85d6SHong Zhang 
951814e85d6SHong Zhang    Logically Collective on TS
952814e85d6SHong Zhang 
953814e85d6SHong Zhang    Input Parameters:
954814e85d6SHong Zhang .  ts - the TS context obtained from TSCreate()
955814e85d6SHong Zhang 
956814e85d6SHong Zhang    Notes:
957814e85d6SHong Zhang    There is no way to remove a single, specific monitor.
958814e85d6SHong Zhang 
959814e85d6SHong Zhang    Level: intermediate
960814e85d6SHong Zhang 
961814e85d6SHong Zhang .keywords: TS, timestep, set, adjoint, monitor
962814e85d6SHong Zhang 
963814e85d6SHong Zhang .seealso: TSAdjointMonitorSet()
964814e85d6SHong Zhang @*/
965814e85d6SHong Zhang PetscErrorCode TSAdjointMonitorCancel(TS ts)
966814e85d6SHong Zhang {
967814e85d6SHong Zhang   PetscErrorCode ierr;
968814e85d6SHong Zhang   PetscInt       i;
969814e85d6SHong Zhang 
970814e85d6SHong Zhang   PetscFunctionBegin;
971814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
972814e85d6SHong Zhang   for (i=0; i<ts->numberadjointmonitors; i++) {
973814e85d6SHong Zhang     if (ts->adjointmonitordestroy[i]) {
974814e85d6SHong Zhang       ierr = (*ts->adjointmonitordestroy[i])(&ts->adjointmonitorcontext[i]);CHKERRQ(ierr);
975814e85d6SHong Zhang     }
976814e85d6SHong Zhang   }
977814e85d6SHong Zhang   ts->numberadjointmonitors = 0;
978814e85d6SHong Zhang   PetscFunctionReturn(0);
979814e85d6SHong Zhang }
980814e85d6SHong Zhang 
981814e85d6SHong Zhang /*@C
982814e85d6SHong Zhang    TSAdjointMonitorDefault - the default monitor of adjoint computations
983814e85d6SHong Zhang 
984814e85d6SHong Zhang    Level: intermediate
985814e85d6SHong Zhang 
986814e85d6SHong Zhang .keywords: TS, set, monitor
987814e85d6SHong Zhang 
988814e85d6SHong Zhang .seealso: TSAdjointMonitorSet()
989814e85d6SHong Zhang @*/
990814e85d6SHong Zhang PetscErrorCode TSAdjointMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,PetscInt numcost,Vec *lambda,Vec *mu,PetscViewerAndFormat *vf)
991814e85d6SHong Zhang {
992814e85d6SHong Zhang   PetscErrorCode ierr;
993814e85d6SHong Zhang   PetscViewer    viewer = vf->viewer;
994814e85d6SHong Zhang 
995814e85d6SHong Zhang   PetscFunctionBegin;
996814e85d6SHong Zhang   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4);
997814e85d6SHong Zhang   ierr = PetscViewerPushFormat(viewer,vf->format);CHKERRQ(ierr);
998814e85d6SHong Zhang   ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
999814e85d6SHong 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);
1000814e85d6SHong Zhang   ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
1001814e85d6SHong Zhang   ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
1002814e85d6SHong Zhang   PetscFunctionReturn(0);
1003814e85d6SHong Zhang }
1004814e85d6SHong Zhang 
1005814e85d6SHong Zhang /*@C
1006814e85d6SHong Zhang    TSAdjointMonitorDrawSensi - Monitors progress of the adjoint TS solvers by calling
1007814e85d6SHong Zhang    VecView() for the sensitivities to initial states at each timestep
1008814e85d6SHong Zhang 
1009814e85d6SHong Zhang    Collective on TS
1010814e85d6SHong Zhang 
1011814e85d6SHong Zhang    Input Parameters:
1012814e85d6SHong Zhang +  ts - the TS context
1013814e85d6SHong Zhang .  step - current time-step
1014814e85d6SHong Zhang .  ptime - current time
1015814e85d6SHong Zhang .  u - current state
1016814e85d6SHong Zhang .  numcost - number of cost functions
1017814e85d6SHong Zhang .  lambda - sensitivities to initial conditions
1018814e85d6SHong Zhang .  mu - sensitivities to parameters
1019814e85d6SHong Zhang -  dummy - either a viewer or NULL
1020814e85d6SHong Zhang 
1021814e85d6SHong Zhang    Level: intermediate
1022814e85d6SHong Zhang 
1023814e85d6SHong Zhang .keywords: TS,  vector, adjoint, monitor, view
1024814e85d6SHong Zhang 
1025814e85d6SHong Zhang .seealso: TSAdjointMonitorSet(), TSAdjointMonitorDefault(), VecView()
1026814e85d6SHong Zhang @*/
1027814e85d6SHong Zhang PetscErrorCode TSAdjointMonitorDrawSensi(TS ts,PetscInt step,PetscReal ptime,Vec u,PetscInt numcost,Vec *lambda,Vec *mu,void *dummy)
1028814e85d6SHong Zhang {
1029814e85d6SHong Zhang   PetscErrorCode   ierr;
1030814e85d6SHong Zhang   TSMonitorDrawCtx ictx = (TSMonitorDrawCtx)dummy;
1031814e85d6SHong Zhang   PetscDraw        draw;
1032814e85d6SHong Zhang   PetscReal        xl,yl,xr,yr,h;
1033814e85d6SHong Zhang   char             time[32];
1034814e85d6SHong Zhang 
1035814e85d6SHong Zhang   PetscFunctionBegin;
1036814e85d6SHong Zhang   if (!(((ictx->howoften > 0) && (!(step % ictx->howoften))) || ((ictx->howoften == -1) && ts->reason))) PetscFunctionReturn(0);
1037814e85d6SHong Zhang 
1038814e85d6SHong Zhang   ierr = VecView(lambda[0],ictx->viewer);CHKERRQ(ierr);
1039814e85d6SHong Zhang   ierr = PetscViewerDrawGetDraw(ictx->viewer,0,&draw);CHKERRQ(ierr);
1040814e85d6SHong Zhang   ierr = PetscSNPrintf(time,32,"Timestep %d Time %g",(int)step,(double)ptime);CHKERRQ(ierr);
1041814e85d6SHong Zhang   ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
1042814e85d6SHong Zhang   h    = yl + .95*(yr - yl);
1043814e85d6SHong Zhang   ierr = PetscDrawStringCentered(draw,.5*(xl+xr),h,PETSC_DRAW_BLACK,time);CHKERRQ(ierr);
1044814e85d6SHong Zhang   ierr = PetscDrawFlush(draw);CHKERRQ(ierr);
1045814e85d6SHong Zhang   PetscFunctionReturn(0);
1046814e85d6SHong Zhang }
1047814e85d6SHong Zhang 
1048814e85d6SHong Zhang /*
1049814e85d6SHong Zhang    TSAdjointSetFromOptions - Sets various TSAdjoint parameters from user options.
1050814e85d6SHong Zhang 
1051814e85d6SHong Zhang    Collective on TSAdjoint
1052814e85d6SHong Zhang 
1053814e85d6SHong Zhang    Input Parameter:
1054814e85d6SHong Zhang .  ts - the TS context
1055814e85d6SHong Zhang 
1056814e85d6SHong Zhang    Options Database Keys:
1057814e85d6SHong Zhang +  -ts_adjoint_solve <yes,no> After solving the ODE/DAE solve the adjoint problem (requires -ts_save_trajectory)
1058814e85d6SHong Zhang .  -ts_adjoint_monitor - print information at each adjoint time step
1059814e85d6SHong Zhang -  -ts_adjoint_monitor_draw_sensi - monitor the sensitivity of the first cost function wrt initial conditions (lambda[0]) graphically
1060814e85d6SHong Zhang 
1061814e85d6SHong Zhang    Level: developer
1062814e85d6SHong Zhang 
1063814e85d6SHong Zhang    Notes:
1064814e85d6SHong Zhang     This is not normally called directly by users
1065814e85d6SHong Zhang 
1066814e85d6SHong Zhang .keywords: TS, trajectory, timestep, set, options, database
1067814e85d6SHong Zhang 
1068814e85d6SHong Zhang .seealso: TSSetSaveTrajectory(), TSTrajectorySetUp()
1069814e85d6SHong Zhang */
1070814e85d6SHong Zhang PetscErrorCode TSAdjointSetFromOptions(PetscOptionItems *PetscOptionsObject,TS ts)
1071814e85d6SHong Zhang {
1072814e85d6SHong Zhang   PetscBool      tflg,opt;
1073814e85d6SHong Zhang   PetscErrorCode ierr;
1074814e85d6SHong Zhang 
1075814e85d6SHong Zhang   PetscFunctionBegin;
1076814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,2);
1077814e85d6SHong Zhang   ierr = PetscOptionsHead(PetscOptionsObject,"TS Adjoint options");CHKERRQ(ierr);
1078814e85d6SHong Zhang   tflg = ts->adjoint_solve ? PETSC_TRUE : PETSC_FALSE;
1079814e85d6SHong Zhang   ierr = PetscOptionsBool("-ts_adjoint_solve","Solve the adjoint problem immediately after solving the forward problem","",tflg,&tflg,&opt);CHKERRQ(ierr);
1080814e85d6SHong Zhang   if (opt) {
1081814e85d6SHong Zhang     ierr = TSSetSaveTrajectory(ts);CHKERRQ(ierr);
1082814e85d6SHong Zhang     ts->adjoint_solve = tflg;
1083814e85d6SHong Zhang   }
1084814e85d6SHong Zhang   ierr = TSAdjointMonitorSetFromOptions(ts,"-ts_adjoint_monitor","Monitor adjoint timestep size","TSAdjointMonitorDefault",TSAdjointMonitorDefault,NULL);CHKERRQ(ierr);
1085814e85d6SHong Zhang   ierr = TSAdjointMonitorSetFromOptions(ts,"-ts_adjoint_monitor_sensi","Monitor sensitivity in the adjoint computation","TSAdjointMonitorSensi",TSAdjointMonitorSensi,NULL);CHKERRQ(ierr);
1086814e85d6SHong Zhang   opt  = PETSC_FALSE;
1087814e85d6SHong Zhang   ierr = PetscOptionsName("-ts_adjoint_monitor_draw_sensi","Monitor adjoint sensitivities (lambda only) graphically","TSAdjointMonitorDrawSensi",&opt);CHKERRQ(ierr);
1088814e85d6SHong Zhang   if (opt) {
1089814e85d6SHong Zhang     TSMonitorDrawCtx ctx;
1090814e85d6SHong Zhang     PetscInt         howoften = 1;
1091814e85d6SHong Zhang 
1092814e85d6SHong Zhang     ierr = PetscOptionsInt("-ts_adjoint_monitor_draw_sensi","Monitor adjoint sensitivities (lambda only) graphically","TSAdjointMonitorDrawSensi",howoften,&howoften,NULL);CHKERRQ(ierr);
1093814e85d6SHong Zhang     ierr = TSMonitorDrawCtxCreate(PetscObjectComm((PetscObject)ts),0,0,PETSC_DECIDE,PETSC_DECIDE,300,300,howoften,&ctx);CHKERRQ(ierr);
1094814e85d6SHong Zhang     ierr = TSAdjointMonitorSet(ts,TSAdjointMonitorDrawSensi,ctx,(PetscErrorCode (*)(void**))TSMonitorDrawCtxDestroy);CHKERRQ(ierr);
1095814e85d6SHong Zhang   }
1096814e85d6SHong Zhang   PetscFunctionReturn(0);
1097814e85d6SHong Zhang }
1098814e85d6SHong Zhang 
1099814e85d6SHong Zhang /*@
1100814e85d6SHong Zhang    TSAdjointStep - Steps one time step backward in the adjoint run
1101814e85d6SHong Zhang 
1102814e85d6SHong Zhang    Collective on TS
1103814e85d6SHong Zhang 
1104814e85d6SHong Zhang    Input Parameter:
1105814e85d6SHong Zhang .  ts - the TS context obtained from TSCreate()
1106814e85d6SHong Zhang 
1107814e85d6SHong Zhang    Level: intermediate
1108814e85d6SHong Zhang 
1109814e85d6SHong Zhang .keywords: TS, adjoint, step
1110814e85d6SHong Zhang 
1111814e85d6SHong Zhang .seealso: TSAdjointSetUp(), TSAdjointSolve()
1112814e85d6SHong Zhang @*/
1113814e85d6SHong Zhang PetscErrorCode TSAdjointStep(TS ts)
1114814e85d6SHong Zhang {
1115814e85d6SHong Zhang   DM               dm;
1116814e85d6SHong Zhang   PetscErrorCode   ierr;
1117814e85d6SHong Zhang 
1118814e85d6SHong Zhang   PetscFunctionBegin;
1119814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1120814e85d6SHong Zhang   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
1121814e85d6SHong Zhang   ierr = TSAdjointSetUp(ts);CHKERRQ(ierr);
1122814e85d6SHong Zhang 
1123814e85d6SHong Zhang   ierr = VecViewFromOptions(ts->vec_sol,(PetscObject)ts,"-ts_view_solution");CHKERRQ(ierr);
1124814e85d6SHong Zhang 
1125814e85d6SHong Zhang   ts->reason = TS_CONVERGED_ITERATING;
1126814e85d6SHong Zhang   ts->ptime_prev = ts->ptime;
1127814e85d6SHong 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);
1128814e85d6SHong Zhang   ierr = PetscLogEventBegin(TS_AdjointStep,ts,0,0,0);CHKERRQ(ierr);
1129814e85d6SHong Zhang   ierr = (*ts->ops->adjointstep)(ts);CHKERRQ(ierr);
1130814e85d6SHong Zhang   ierr = PetscLogEventEnd(TS_AdjointStep,ts,0,0,0);CHKERRQ(ierr);
1131814e85d6SHong Zhang   ts->adjoint_steps++; ts->steps--;
1132814e85d6SHong Zhang 
1133814e85d6SHong Zhang   if (ts->reason < 0) {
1134814e85d6SHong Zhang     if (ts->errorifstepfailed) {
1135814e85d6SHong Zhang       if (ts->reason == TS_DIVERGED_NONLINEAR_SOLVE) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_NOT_CONVERGED,"TSStep has failed due to %s, increase -ts_max_snes_failures or make negative to attempt recovery",TSConvergedReasons[ts->reason]);
1136814e85d6SHong Zhang       else if (ts->reason == TS_DIVERGED_STEP_REJECTED) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_NOT_CONVERGED,"TSStep has failed due to %s, increase -ts_max_reject or make negative to attempt recovery",TSConvergedReasons[ts->reason]);
1137814e85d6SHong Zhang       else SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_NOT_CONVERGED,"TSStep has failed due to %s",TSConvergedReasons[ts->reason]);
1138814e85d6SHong Zhang     }
1139814e85d6SHong Zhang   } else if (!ts->reason) {
1140814e85d6SHong Zhang     if (ts->adjoint_steps >= ts->adjoint_max_steps) ts->reason = TS_CONVERGED_ITS;
1141814e85d6SHong Zhang   }
1142814e85d6SHong Zhang   PetscFunctionReturn(0);
1143814e85d6SHong Zhang }
1144814e85d6SHong Zhang 
1145814e85d6SHong Zhang /*@
1146814e85d6SHong Zhang    TSAdjointSolve - Solves the discrete ajoint problem for an ODE/DAE
1147814e85d6SHong Zhang 
1148814e85d6SHong Zhang    Collective on TS
1149814e85d6SHong Zhang 
1150814e85d6SHong Zhang    Input Parameter:
1151814e85d6SHong Zhang .  ts - the TS context obtained from TSCreate()
1152814e85d6SHong Zhang 
1153814e85d6SHong Zhang    Options Database:
1154814e85d6SHong Zhang . -ts_adjoint_view_solution <viewerinfo> - views the first gradient with respect to the initial values
1155814e85d6SHong Zhang 
1156814e85d6SHong Zhang    Level: intermediate
1157814e85d6SHong Zhang 
1158814e85d6SHong Zhang    Notes:
1159814e85d6SHong Zhang    This must be called after a call to TSSolve() that solves the forward problem
1160814e85d6SHong Zhang 
1161814e85d6SHong Zhang    By default this will integrate back to the initial time, one can use TSAdjointSetSteps() to step back to a later time
1162814e85d6SHong Zhang 
1163814e85d6SHong Zhang .keywords: TS, timestep, solve
1164814e85d6SHong Zhang 
1165814e85d6SHong Zhang .seealso: TSCreate(), TSSetCostGradients(), TSSetSolution(), TSAdjointStep()
1166814e85d6SHong Zhang @*/
1167814e85d6SHong Zhang PetscErrorCode TSAdjointSolve(TS ts)
1168814e85d6SHong Zhang {
1169814e85d6SHong Zhang   PetscErrorCode    ierr;
1170814e85d6SHong Zhang 
1171814e85d6SHong Zhang   PetscFunctionBegin;
1172814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1173814e85d6SHong Zhang   ierr = TSAdjointSetUp(ts);CHKERRQ(ierr);
1174814e85d6SHong Zhang 
1175814e85d6SHong Zhang   /* reset time step and iteration counters */
1176814e85d6SHong Zhang   ts->adjoint_steps     = 0;
1177814e85d6SHong Zhang   ts->ksp_its           = 0;
1178814e85d6SHong Zhang   ts->snes_its          = 0;
1179814e85d6SHong Zhang   ts->num_snes_failures = 0;
1180814e85d6SHong Zhang   ts->reject            = 0;
1181814e85d6SHong Zhang   ts->reason            = TS_CONVERGED_ITERATING;
1182814e85d6SHong Zhang 
1183814e85d6SHong Zhang   if (!ts->adjoint_max_steps) ts->adjoint_max_steps = ts->steps;
1184814e85d6SHong Zhang   if (ts->adjoint_steps >= ts->adjoint_max_steps) ts->reason = TS_CONVERGED_ITS;
1185814e85d6SHong Zhang 
1186814e85d6SHong Zhang   while (!ts->reason) {
1187814e85d6SHong Zhang     ierr = TSTrajectoryGet(ts->trajectory,ts,ts->steps,&ts->ptime);CHKERRQ(ierr);
1188814e85d6SHong Zhang     ierr = TSAdjointMonitor(ts,ts->steps,ts->ptime,ts->vec_sol,ts->numcost,ts->vecs_sensi,ts->vecs_sensip);CHKERRQ(ierr);
1189814e85d6SHong Zhang     ierr = TSAdjointEventHandler(ts);CHKERRQ(ierr);
1190814e85d6SHong Zhang     ierr = TSAdjointStep(ts);CHKERRQ(ierr);
1191814e85d6SHong Zhang     if (ts->vec_costintegral && !ts->costintegralfwd) {
1192814e85d6SHong Zhang       ierr = TSAdjointCostIntegral(ts);CHKERRQ(ierr);
1193814e85d6SHong Zhang     }
1194814e85d6SHong Zhang   }
1195814e85d6SHong Zhang   ierr = TSTrajectoryGet(ts->trajectory,ts,ts->steps,&ts->ptime);CHKERRQ(ierr);
1196814e85d6SHong Zhang   ierr = TSAdjointMonitor(ts,ts->steps,ts->ptime,ts->vec_sol,ts->numcost,ts->vecs_sensi,ts->vecs_sensip);CHKERRQ(ierr);
1197814e85d6SHong Zhang   ts->solvetime = ts->ptime;
1198814e85d6SHong Zhang   ierr = TSTrajectoryViewFromOptions(ts->trajectory,NULL,"-ts_trajectory_view");CHKERRQ(ierr);
1199814e85d6SHong Zhang   ierr = VecViewFromOptions(ts->vecs_sensi[0],(PetscObject) ts, "-ts_adjoint_view_solution");CHKERRQ(ierr);
1200814e85d6SHong Zhang   ts->adjoint_max_steps = 0;
1201814e85d6SHong Zhang   PetscFunctionReturn(0);
1202814e85d6SHong Zhang }
1203814e85d6SHong Zhang 
1204814e85d6SHong Zhang /*@C
1205814e85d6SHong Zhang    TSAdjointMonitor - Runs all user-provided adjoint monitor routines set using TSAdjointMonitorSet()
1206814e85d6SHong Zhang 
1207814e85d6SHong Zhang    Collective on TS
1208814e85d6SHong Zhang 
1209814e85d6SHong Zhang    Input Parameters:
1210814e85d6SHong Zhang +  ts - time stepping context obtained from TSCreate()
1211814e85d6SHong Zhang .  step - step number that has just completed
1212814e85d6SHong Zhang .  ptime - model time of the state
1213814e85d6SHong Zhang .  u - state at the current model time
1214814e85d6SHong Zhang .  numcost - number of cost functions (dimension of lambda  or mu)
1215814e85d6SHong Zhang .  lambda - vectors containing the gradients of the cost functions with respect to the ODE/DAE solution variables
1216814e85d6SHong Zhang -  mu - vectors containing the gradients of the cost functions with respect to the problem parameters
1217814e85d6SHong Zhang 
1218814e85d6SHong Zhang    Notes:
1219814e85d6SHong Zhang    TSAdjointMonitor() is typically used automatically within the time stepping implementations.
1220814e85d6SHong Zhang    Users would almost never call this routine directly.
1221814e85d6SHong Zhang 
1222814e85d6SHong Zhang    Level: developer
1223814e85d6SHong Zhang 
1224814e85d6SHong Zhang .keywords: TS, timestep
1225814e85d6SHong Zhang @*/
1226814e85d6SHong Zhang PetscErrorCode TSAdjointMonitor(TS ts,PetscInt step,PetscReal ptime,Vec u,PetscInt numcost,Vec *lambda, Vec *mu)
1227814e85d6SHong Zhang {
1228814e85d6SHong Zhang   PetscErrorCode ierr;
1229814e85d6SHong Zhang   PetscInt       i,n = ts->numberadjointmonitors;
1230814e85d6SHong Zhang 
1231814e85d6SHong Zhang   PetscFunctionBegin;
1232814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1233814e85d6SHong Zhang   PetscValidHeaderSpecific(u,VEC_CLASSID,4);
12348860a134SJunchao Zhang   ierr = VecLockReadPush(u);CHKERRQ(ierr);
1235814e85d6SHong Zhang   for (i=0; i<n; i++) {
1236814e85d6SHong Zhang     ierr = (*ts->adjointmonitor[i])(ts,step,ptime,u,numcost,lambda,mu,ts->adjointmonitorcontext[i]);CHKERRQ(ierr);
1237814e85d6SHong Zhang   }
12388860a134SJunchao Zhang   ierr = VecLockReadPop(u);CHKERRQ(ierr);
1239814e85d6SHong Zhang   PetscFunctionReturn(0);
1240814e85d6SHong Zhang }
1241814e85d6SHong Zhang 
1242814e85d6SHong Zhang /*@
1243814e85d6SHong Zhang  TSAdjointCostIntegral - Evaluate the cost integral in the adjoint run.
1244814e85d6SHong Zhang 
1245814e85d6SHong Zhang  Collective on TS
1246814e85d6SHong Zhang 
1247814e85d6SHong Zhang  Input Arguments:
1248814e85d6SHong Zhang  .  ts - time stepping context
1249814e85d6SHong Zhang 
1250814e85d6SHong Zhang  Level: advanced
1251814e85d6SHong Zhang 
1252814e85d6SHong Zhang  Notes:
1253814e85d6SHong Zhang  This function cannot be called until TSAdjointStep() has been completed.
1254814e85d6SHong Zhang 
1255814e85d6SHong Zhang  .seealso: TSAdjointSolve(), TSAdjointStep
1256814e85d6SHong Zhang  @*/
1257814e85d6SHong Zhang PetscErrorCode TSAdjointCostIntegral(TS ts)
1258814e85d6SHong Zhang {
1259814e85d6SHong Zhang     PetscErrorCode ierr;
1260814e85d6SHong Zhang     PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1261814e85d6SHong 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);
1262814e85d6SHong Zhang     ierr = (*ts->ops->adjointintegral)(ts);CHKERRQ(ierr);
1263814e85d6SHong Zhang     PetscFunctionReturn(0);
1264814e85d6SHong Zhang }
1265814e85d6SHong Zhang 
1266814e85d6SHong Zhang /* ------------------ Forward (tangent linear) sensitivity  ------------------*/
1267814e85d6SHong Zhang 
1268814e85d6SHong Zhang /*@
1269814e85d6SHong Zhang   TSForwardSetUp - Sets up the internal data structures for the later use
1270814e85d6SHong Zhang   of forward sensitivity analysis
1271814e85d6SHong Zhang 
1272814e85d6SHong Zhang   Collective on TS
1273814e85d6SHong Zhang 
1274814e85d6SHong Zhang   Input Parameter:
1275814e85d6SHong Zhang . ts - the TS context obtained from TSCreate()
1276814e85d6SHong Zhang 
1277814e85d6SHong Zhang   Level: advanced
1278814e85d6SHong Zhang 
1279814e85d6SHong Zhang .keywords: TS, forward sensitivity, setup
1280814e85d6SHong Zhang 
1281814e85d6SHong Zhang .seealso: TSCreate(), TSDestroy(), TSSetUp()
1282814e85d6SHong Zhang @*/
1283814e85d6SHong Zhang PetscErrorCode TSForwardSetUp(TS ts)
1284814e85d6SHong Zhang {
1285814e85d6SHong Zhang   PetscErrorCode ierr;
1286814e85d6SHong Zhang 
1287814e85d6SHong Zhang   PetscFunctionBegin;
1288814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1289814e85d6SHong Zhang   if (ts->forwardsetupcalled) PetscFunctionReturn(0);
1290814e85d6SHong Zhang   if (ts->vec_costintegral && !ts->vecs_integral_sensip ) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TSForwardSetIntegralGradients() before TSSetCostIntegrand()");
1291814e85d6SHong Zhang   if (ts->vecs_integral_sensip) {
1292c9ad9de0SHong Zhang     ierr = VecDuplicateVecs(ts->vec_sol,ts->numcost,&ts->vecs_drdu);CHKERRQ(ierr);
1293814e85d6SHong Zhang     ierr = VecDuplicateVecs(ts->vecs_integral_sensip[0],ts->numcost,&ts->vecs_drdp);CHKERRQ(ierr);
1294814e85d6SHong Zhang   }
1295814e85d6SHong Zhang 
1296814e85d6SHong Zhang   if (ts->ops->forwardsetup) {
1297814e85d6SHong Zhang     ierr = (*ts->ops->forwardsetup)(ts);CHKERRQ(ierr);
1298814e85d6SHong Zhang   }
1299814e85d6SHong Zhang   ts->forwardsetupcalled = PETSC_TRUE;
1300814e85d6SHong Zhang   PetscFunctionReturn(0);
1301814e85d6SHong Zhang }
1302814e85d6SHong Zhang 
1303814e85d6SHong Zhang /*@
1304*7adebddeSHong Zhang   TSForwardReset - Reset the internal data structures used by forward sensitivity analysis
1305*7adebddeSHong Zhang 
1306*7adebddeSHong Zhang   Collective on TS
1307*7adebddeSHong Zhang 
1308*7adebddeSHong Zhang   Input Parameter:
1309*7adebddeSHong Zhang . ts - the TS context obtained from TSCreate()
1310*7adebddeSHong Zhang 
1311*7adebddeSHong Zhang   Level: advanced
1312*7adebddeSHong Zhang 
1313*7adebddeSHong Zhang .keywords: TS, forward sensitivity, reset
1314*7adebddeSHong Zhang 
1315*7adebddeSHong Zhang .seealso: TSCreate(), TSDestroy(), TSForwardSetUp()
1316*7adebddeSHong Zhang @*/
1317*7adebddeSHong Zhang PetscErrorCode TSForwardReset(TS ts)
1318*7adebddeSHong Zhang {
1319*7adebddeSHong Zhang   PetscErrorCode ierr;
1320*7adebddeSHong Zhang 
1321*7adebddeSHong Zhang   PetscFunctionBegin;
1322*7adebddeSHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1323*7adebddeSHong Zhang   if (ts->ops->forwardreset) {
1324*7adebddeSHong Zhang     ierr = (*ts->ops->forwardreset)(ts);CHKERRQ(ierr);
1325*7adebddeSHong Zhang   }
1326*7adebddeSHong Zhang   ierr = MatDestroy(&ts->mat_sensip);CHKERRQ(ierr);
1327*7adebddeSHong Zhang   ts->vecs_integral_sensip = NULL;
1328*7adebddeSHong Zhang   ts->forward_solve        = PETSC_FALSE;
1329*7adebddeSHong Zhang   ts->forwardsetupcalled   = PETSC_FALSE;
1330*7adebddeSHong Zhang   PetscFunctionReturn(0);
1331*7adebddeSHong Zhang }
1332*7adebddeSHong Zhang 
1333*7adebddeSHong Zhang /*@
1334814e85d6SHong Zhang   TSForwardSetIntegralGradients - Set the vectors holding forward sensitivities of the integral term.
1335814e85d6SHong Zhang 
1336814e85d6SHong Zhang   Input Parameter:
1337814e85d6SHong Zhang . ts- the TS context obtained from TSCreate()
1338814e85d6SHong Zhang . numfwdint- number of integrals
1339814e85d6SHong Zhang . vp = the vectors containing the gradients for each integral w.r.t. parameters
1340814e85d6SHong Zhang 
1341814e85d6SHong Zhang   Level: intermediate
1342814e85d6SHong Zhang 
1343814e85d6SHong Zhang .keywords: TS, forward sensitivity
1344814e85d6SHong Zhang 
1345814e85d6SHong Zhang .seealso: TSForwardGetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep()
1346814e85d6SHong Zhang @*/
1347814e85d6SHong Zhang PetscErrorCode TSForwardSetIntegralGradients(TS ts,PetscInt numfwdint,Vec *vp)
1348814e85d6SHong Zhang {
1349814e85d6SHong Zhang   PetscFunctionBegin;
1350814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1351814e85d6SHong 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()");
1352814e85d6SHong Zhang   if (!ts->numcost) ts->numcost = numfwdint;
1353814e85d6SHong Zhang 
1354814e85d6SHong Zhang   ts->vecs_integral_sensip = vp;
1355814e85d6SHong Zhang   PetscFunctionReturn(0);
1356814e85d6SHong Zhang }
1357814e85d6SHong Zhang 
1358814e85d6SHong Zhang /*@
1359814e85d6SHong Zhang   TSForwardGetIntegralGradients - Returns the forward sensitivities ofthe integral term.
1360814e85d6SHong Zhang 
1361814e85d6SHong Zhang   Input Parameter:
1362814e85d6SHong Zhang . ts- the TS context obtained from TSCreate()
1363814e85d6SHong Zhang 
1364814e85d6SHong Zhang   Output Parameter:
1365814e85d6SHong Zhang . vp = the vectors containing the gradients for each integral w.r.t. parameters
1366814e85d6SHong Zhang 
1367814e85d6SHong Zhang   Level: intermediate
1368814e85d6SHong Zhang 
1369814e85d6SHong Zhang .keywords: TS, forward sensitivity
1370814e85d6SHong Zhang 
1371814e85d6SHong Zhang .seealso: TSForwardSetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep()
1372814e85d6SHong Zhang @*/
1373814e85d6SHong Zhang PetscErrorCode TSForwardGetIntegralGradients(TS ts,PetscInt *numfwdint,Vec **vp)
1374814e85d6SHong Zhang {
1375814e85d6SHong Zhang   PetscFunctionBegin;
1376814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1377814e85d6SHong Zhang   PetscValidPointer(vp,3);
1378814e85d6SHong Zhang   if (numfwdint) *numfwdint = ts->numcost;
1379814e85d6SHong Zhang   if (vp) *vp = ts->vecs_integral_sensip;
1380814e85d6SHong Zhang   PetscFunctionReturn(0);
1381814e85d6SHong Zhang }
1382814e85d6SHong Zhang 
1383814e85d6SHong Zhang /*@
1384814e85d6SHong Zhang   TSForwardStep - Compute the forward sensitivity for one time step.
1385814e85d6SHong Zhang 
1386814e85d6SHong Zhang   Collective on TS
1387814e85d6SHong Zhang 
1388814e85d6SHong Zhang   Input Arguments:
1389814e85d6SHong Zhang . ts - time stepping context
1390814e85d6SHong Zhang 
1391814e85d6SHong Zhang   Level: advanced
1392814e85d6SHong Zhang 
1393814e85d6SHong Zhang   Notes:
1394814e85d6SHong Zhang   This function cannot be called until TSStep() has been completed.
1395814e85d6SHong Zhang 
1396814e85d6SHong Zhang .keywords: TS, forward sensitivity
1397814e85d6SHong Zhang 
1398814e85d6SHong Zhang .seealso: TSForwardSetSensitivities(), TSForwardGetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardSetUp()
1399814e85d6SHong Zhang @*/
1400814e85d6SHong Zhang PetscErrorCode TSForwardStep(TS ts)
1401814e85d6SHong Zhang {
1402814e85d6SHong Zhang   PetscErrorCode ierr;
1403814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1404814e85d6SHong Zhang   if (!ts->ops->forwardstep) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"%s does not provide forward sensitivity analysis",((PetscObject)ts)->type_name);
1405814e85d6SHong Zhang   ierr = PetscLogEventBegin(TS_ForwardStep,ts,0,0,0);CHKERRQ(ierr);
1406814e85d6SHong Zhang   ierr = (*ts->ops->forwardstep)(ts);CHKERRQ(ierr);
1407814e85d6SHong Zhang   ierr = PetscLogEventEnd(TS_ForwardStep,ts,0,0,0);CHKERRQ(ierr);
1408814e85d6SHong Zhang   PetscFunctionReturn(0);
1409814e85d6SHong Zhang }
1410814e85d6SHong Zhang 
1411814e85d6SHong Zhang /*@
1412814e85d6SHong Zhang   TSForwardSetSensitivities - Sets the initial value of the trajectory sensitivities of solution  w.r.t. the problem parameters and initial values.
1413814e85d6SHong Zhang 
1414814e85d6SHong Zhang   Logically Collective on TS and Vec
1415814e85d6SHong Zhang 
1416814e85d6SHong Zhang   Input Parameters:
1417814e85d6SHong Zhang + ts - the TS context obtained from TSCreate()
1418814e85d6SHong Zhang . nump - number of parameters
1419814e85d6SHong Zhang - Smat - sensitivities with respect to the parameters, the number of entries in these vectors is the same as the number of parameters
1420814e85d6SHong Zhang 
1421814e85d6SHong Zhang   Level: beginner
1422814e85d6SHong Zhang 
1423814e85d6SHong Zhang   Notes:
1424814e85d6SHong Zhang   Forward sensitivity is also called 'trajectory sensitivity' in some fields such as power systems.
1425814e85d6SHong Zhang   This function turns on a flag to trigger TSSolve() to compute forward sensitivities automatically.
1426814e85d6SHong Zhang   You must call this function before TSSolve().
1427814e85d6SHong Zhang   The entries in the sensitivity matrix must be correctly initialized with the values S = dy/dp|startingtime.
1428814e85d6SHong Zhang 
1429814e85d6SHong Zhang .keywords: TS, timestep, set, forward sensitivity, initial values
1430814e85d6SHong Zhang 
1431814e85d6SHong Zhang .seealso: TSForwardGetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep()
1432814e85d6SHong Zhang @*/
1433814e85d6SHong Zhang PetscErrorCode TSForwardSetSensitivities(TS ts,PetscInt nump,Mat Smat)
1434814e85d6SHong Zhang {
1435814e85d6SHong Zhang   PetscErrorCode ierr;
1436814e85d6SHong Zhang 
1437814e85d6SHong Zhang   PetscFunctionBegin;
1438814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1439814e85d6SHong Zhang   PetscValidHeaderSpecific(Smat,MAT_CLASSID,3);
1440814e85d6SHong Zhang   ts->forward_solve  = PETSC_TRUE;
1441b5b4017aSHong Zhang   if (nump == PETSC_DEFAULT) {
1442b5b4017aSHong Zhang     ierr = MatGetSize(Smat,NULL,&ts->num_parameters);CHKERRQ(ierr);
1443b5b4017aSHong Zhang   } else ts->num_parameters = nump;
1444814e85d6SHong Zhang   ierr = PetscObjectReference((PetscObject)Smat);CHKERRQ(ierr);
1445814e85d6SHong Zhang   ierr = MatDestroy(&ts->mat_sensip);CHKERRQ(ierr);
1446814e85d6SHong Zhang   ts->mat_sensip = Smat;
1447814e85d6SHong Zhang   PetscFunctionReturn(0);
1448814e85d6SHong Zhang }
1449814e85d6SHong Zhang 
1450814e85d6SHong Zhang /*@
1451814e85d6SHong Zhang   TSForwardGetSensitivities - Returns the trajectory sensitivities
1452814e85d6SHong Zhang 
1453814e85d6SHong Zhang   Not Collective, but Vec returned is parallel if TS is parallel
1454814e85d6SHong Zhang 
1455814e85d6SHong Zhang   Output Parameter:
1456814e85d6SHong Zhang + ts - the TS context obtained from TSCreate()
1457814e85d6SHong Zhang . nump - number of parameters
1458814e85d6SHong Zhang - Smat - sensitivities with respect to the parameters, the number of entries in these vectors is the same as the number of parameters
1459814e85d6SHong Zhang 
1460814e85d6SHong Zhang   Level: intermediate
1461814e85d6SHong Zhang 
1462814e85d6SHong Zhang .keywords: TS, forward sensitivity
1463814e85d6SHong Zhang 
1464814e85d6SHong Zhang .seealso: TSForwardSetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep()
1465814e85d6SHong Zhang @*/
1466814e85d6SHong Zhang PetscErrorCode TSForwardGetSensitivities(TS ts,PetscInt *nump,Mat *Smat)
1467814e85d6SHong Zhang {
1468814e85d6SHong Zhang   PetscFunctionBegin;
1469814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1470814e85d6SHong Zhang   if (nump) *nump = ts->num_parameters;
1471814e85d6SHong Zhang   if (Smat) *Smat = ts->mat_sensip;
1472814e85d6SHong Zhang   PetscFunctionReturn(0);
1473814e85d6SHong Zhang }
1474814e85d6SHong Zhang 
1475814e85d6SHong Zhang /*@
1476814e85d6SHong Zhang    TSForwardCostIntegral - Evaluate the cost integral in the forward run.
1477814e85d6SHong Zhang 
1478814e85d6SHong Zhang    Collective on TS
1479814e85d6SHong Zhang 
1480814e85d6SHong Zhang    Input Arguments:
1481814e85d6SHong Zhang .  ts - time stepping context
1482814e85d6SHong Zhang 
1483814e85d6SHong Zhang    Level: advanced
1484814e85d6SHong Zhang 
1485814e85d6SHong Zhang    Notes:
1486814e85d6SHong Zhang    This function cannot be called until TSStep() has been completed.
1487814e85d6SHong Zhang 
1488814e85d6SHong Zhang .seealso: TSSolve(), TSAdjointCostIntegral()
1489814e85d6SHong Zhang @*/
1490814e85d6SHong Zhang PetscErrorCode TSForwardCostIntegral(TS ts)
1491814e85d6SHong Zhang {
1492814e85d6SHong Zhang   PetscErrorCode ierr;
1493814e85d6SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1494814e85d6SHong 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);
1495814e85d6SHong Zhang   ierr = (*ts->ops->forwardintegral)(ts);CHKERRQ(ierr);
1496814e85d6SHong Zhang   PetscFunctionReturn(0);
1497814e85d6SHong Zhang }
1498b5b4017aSHong Zhang 
1499b5b4017aSHong Zhang /*@
1500b5b4017aSHong Zhang   TSForwardSetInitialSensitivities - Set initial values for tangent linear sensitivities
1501b5b4017aSHong Zhang 
1502b5b4017aSHong Zhang   Collective on TS and Mat
1503b5b4017aSHong Zhang 
1504b5b4017aSHong Zhang   Input Parameter
1505b5b4017aSHong Zhang + ts - the TS context obtained from TSCreate()
1506b5b4017aSHong Zhang - didp - parametric sensitivities of the initial condition
1507b5b4017aSHong Zhang 
1508b5b4017aSHong Zhang   Level: intermediate
1509b5b4017aSHong Zhang 
1510b5b4017aSHong 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.
1511b5b4017aSHong Zhang 
1512b5b4017aSHong Zhang .seealso: TSForwardSetSensitivities()
1513b5b4017aSHong Zhang @*/
1514b5b4017aSHong Zhang PetscErrorCode TSForwardSetInitialSensitivities(TS ts,Mat didp)
1515b5b4017aSHong Zhang {
1516b5b4017aSHong Zhang   Vec            sp;
1517b5b4017aSHong Zhang   PetscInt       lsize;
1518b5b4017aSHong Zhang   PetscScalar    *xarr;
1519b5b4017aSHong Zhang   PetscErrorCode ierr;
1520b5b4017aSHong Zhang 
1521b5b4017aSHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1522b5b4017aSHong Zhang   if (ts->vec_dir) { /* indicates second-order adjoint caculation */
15236affb6f8SHong Zhang     Mat A;
15246affb6f8SHong Zhang     ierr = TSForwardGetSensitivities(ts,NULL,&A);CHKERRQ(ierr);
15256affb6f8SHong Zhang     if (!A) { /* create a single-column dense matrix */
1526b5b4017aSHong Zhang       ierr = VecGetLocalSize(ts->vec_dir,&lsize);CHKERRQ(ierr);
15276affb6f8SHong Zhang       ierr = MatCreateDense(PETSC_COMM_WORLD,lsize,PETSC_DECIDE,PETSC_DECIDE,1,NULL,&A);CHKERRQ(ierr);
1528b5b4017aSHong Zhang     }
1529b5b4017aSHong Zhang     ierr = VecDuplicate(ts->vec_dir,&sp);CHKERRQ(ierr);
15306affb6f8SHong Zhang     ierr = MatDenseGetColumn(A,0,&xarr);CHKERRQ(ierr);
1531b5b4017aSHong Zhang     ierr = VecPlaceArray(sp,xarr);CHKERRQ(ierr);
1532b5b4017aSHong Zhang     if (didp) {
1533b5b4017aSHong Zhang       ierr = MatMult(didp,ts->vec_dir,sp);CHKERRQ(ierr);
1534b5b4017aSHong Zhang     } else { /* identity matrix assumed */
1535b5b4017aSHong Zhang       ierr = VecCopy(ts->vec_dir,sp);CHKERRQ(ierr);
1536b5b4017aSHong Zhang     }
1537b5b4017aSHong Zhang     ierr = VecResetArray(sp);CHKERRQ(ierr);
15386affb6f8SHong Zhang     ierr = MatDenseRestoreColumn(A,&xarr);CHKERRQ(ierr);
1539b5b4017aSHong Zhang     ierr = VecDestroy(&sp);CHKERRQ(ierr);
15406affb6f8SHong Zhang     ierr = TSForwardSetSensitivities(ts,1,A);CHKERRQ(ierr);
1541b5b4017aSHong Zhang   } else {
1542b5b4017aSHong Zhang     PetscValidHeaderSpecific(didp,MAT_CLASSID,2);
1543b5b4017aSHong Zhang     if (!ts->mat_sensip) {
1544b5b4017aSHong Zhang       ierr = TSForwardSetSensitivities(ts,PETSC_DEFAULT,didp);CHKERRQ(ierr);
1545b5b4017aSHong Zhang     }
1546b5b4017aSHong Zhang   }
1547b5b4017aSHong Zhang   PetscFunctionReturn(0);
1548b5b4017aSHong Zhang }
15496affb6f8SHong Zhang 
15506affb6f8SHong Zhang /*@
15516affb6f8SHong Zhang    TSForwardGetStages - Get the number of stages and the tangent linear sensitivities at the intermediate stages
15526affb6f8SHong Zhang 
15536affb6f8SHong Zhang    Input Parameter:
15546affb6f8SHong Zhang .  ts - the TS context obtained from TSCreate()
15556affb6f8SHong Zhang 
15566affb6f8SHong Zhang    Output Parameters:
15576affb6f8SHong Zhang +  ns - nu
15586affb6f8SHong Zhang -  S - tangent linear sensitivities at the intermediate stages
15596affb6f8SHong Zhang 
15606affb6f8SHong Zhang    Level: advanced
15616affb6f8SHong Zhang 
15626affb6f8SHong Zhang .keywords: TS, second-order adjoint, forward sensitivity
15636affb6f8SHong Zhang @*/
15646affb6f8SHong Zhang PetscErrorCode TSForwardGetStages(TS ts,PetscInt *ns,Mat **S)
15656affb6f8SHong Zhang {
15666affb6f8SHong Zhang   PetscErrorCode ierr;
15676affb6f8SHong Zhang 
15686affb6f8SHong Zhang   PetscFunctionBegin;
15696affb6f8SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
15706affb6f8SHong Zhang 
15716affb6f8SHong Zhang   if (!ts->ops->getstages) *S=NULL;
15726affb6f8SHong Zhang   else {
15736affb6f8SHong Zhang     ierr = (*ts->ops->forwardgetstages)(ts,ns,S);CHKERRQ(ierr);
15746affb6f8SHong Zhang   }
15756affb6f8SHong Zhang   PetscFunctionReturn(0);
15766affb6f8SHong Zhang }
1577