Actual source code: tssen.c
petsc-3.12.5 2020-03-29
1: #include <petsc/private/tsimpl.h>
2: #include <petscdraw.h>
4: PetscLogEvent TS_AdjointStep,TS_ForwardStep,TS_JacobianPEval;
6: /* #define TSADJOINT_STAGE */
8: /* ------------------------ Sensitivity Context ---------------------------*/
10: /*@C
11: 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.
13: Logically Collective on TS
15: Input Parameters:
16: + ts - TS context obtained from TSCreate()
17: . Amat - JacobianP matrix
18: . func - function
19: - ctx - [optional] user-defined function context
21: Calling sequence of func:
22: $ func (TS ts,PetscReal t,Vec y,Mat A,void *ctx);
23: + t - current timestep
24: . U - input vector (current ODE solution)
25: . A - output matrix
26: - ctx - [optional] user-defined function context
28: Level: intermediate
30: Notes:
31: 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
33: .seealso: TSGetRHSJacobianP()
34: @*/
35: PetscErrorCode TSSetRHSJacobianP(TS ts,Mat Amat,PetscErrorCode (*func)(TS,PetscReal,Vec,Mat,void*),void *ctx)
36: {
43: ts->rhsjacobianp = func;
44: ts->rhsjacobianpctx = ctx;
45: if(Amat) {
46: PetscObjectReference((PetscObject)Amat);
47: MatDestroy(&ts->Jacprhs);
48: ts->Jacprhs = Amat;
49: }
50: return(0);
51: }
53: /*@C
54: TSGetRHSJacobianP - Gets the function that computes the Jacobian of G w.r.t. the parameters P where U_t = G(U,P,t), as well as the location to store the matrix.
56: Logically Collective on TS
58: Input Parameters:
59: . ts - TS context obtained from TSCreate()
61: Output Parameters:
62: + Amat - JacobianP matrix
63: . func - function
64: - ctx - [optional] user-defined function context
66: Calling sequence of func:
67: $ func (TS ts,PetscReal t,Vec y,Mat A,void *ctx);
68: + t - current timestep
69: . U - input vector (current ODE solution)
70: . A - output matrix
71: - ctx - [optional] user-defined function context
73: Level: intermediate
75: Notes:
76: 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
78: .seealso: TSSetRHSJacobianP()
79: @*/
80: PetscErrorCode TSGetRHSJacobianP(TS ts,Mat *Amat,PetscErrorCode (**func)(TS,PetscReal,Vec,Mat,void*),void **ctx)
81: {
83: if (func) *func = ts->rhsjacobianp;
84: if (ctx) *ctx = ts->rhsjacobianpctx;
85: if (Amat) *Amat = ts->Jacprhs;
86: return(0);
87: }
89: /*@C
90: TSComputeRHSJacobianP - Runs the user-defined JacobianP function.
92: Collective on TS
94: Input Parameters:
95: . ts - The TS context obtained from TSCreate()
97: Level: developer
99: .seealso: TSSetRHSJacobianP()
100: @*/
101: PetscErrorCode TSComputeRHSJacobianP(TS ts,PetscReal t,Vec U,Mat Amat)
102: {
106: if (!Amat) return(0);
110: PetscStackPush("TS user JacobianP function for sensitivity analysis");
111: (*ts->rhsjacobianp)(ts,t,U,Amat,ts->rhsjacobianpctx);
112: PetscStackPop;
113: return(0);
114: }
116: /*@C
117: 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.
119: Logically Collective on TS
121: Input Parameters:
122: + ts - TS context obtained from TSCreate()
123: . Amat - JacobianP matrix
124: . func - function
125: - ctx - [optional] user-defined function context
127: Calling sequence of func:
128: $ func (TS ts,PetscReal t,Vec y,Mat A,void *ctx);
129: + t - current timestep
130: . U - input vector (current ODE solution)
131: . Udot - time derivative of state vector
132: . shift - shift to apply, see note below
133: . A - output matrix
134: - ctx - [optional] user-defined function context
136: Level: intermediate
138: Notes:
139: 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
141: .seealso:
142: @*/
143: PetscErrorCode TSSetIJacobianP(TS ts,Mat Amat,PetscErrorCode (*func)(TS,PetscReal,Vec,Vec,PetscReal,Mat,void*),void *ctx)
144: {
151: ts->ijacobianp = func;
152: ts->ijacobianpctx = ctx;
153: if(Amat) {
154: PetscObjectReference((PetscObject)Amat);
155: MatDestroy(&ts->Jacp);
156: ts->Jacp = Amat;
157: }
158: return(0);
159: }
161: /*@C
162: TSComputeIJacobianP - Runs the user-defined IJacobianP function.
164: Collective on TS
166: Input Parameters:
167: + ts - the TS context
168: . t - current timestep
169: . U - state vector
170: . Udot - time derivative of state vector
171: . shift - shift to apply, see note below
172: - imex - flag indicates if the method is IMEX so that the RHSJacobian should be kept separate
174: Output Parameters:
175: . A - Jacobian matrix
177: Level: developer
179: .seealso: TSSetIJacobianP()
180: @*/
181: PetscErrorCode TSComputeIJacobianP(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal shift,Mat Amat,PetscBool imex)
182: {
186: if (!Amat) return(0);
191: PetscLogEventBegin(TS_JacobianPEval,ts,U,Amat,0);
192: if (ts->ijacobianp) {
193: PetscStackPush("TS user JacobianP function for sensitivity analysis");
194: (*ts->ijacobianp)(ts,t,U,Udot,shift,Amat,ts->ijacobianpctx);
195: PetscStackPop;
196: }
197: if (imex) {
198: if (!ts->ijacobianp) { /* system was written as Udot = G(t,U) */
199: PetscBool assembled;
200: MatZeroEntries(Amat);
201: MatAssembled(Amat,&assembled);
202: if (!assembled) {
203: MatAssemblyBegin(Amat,MAT_FINAL_ASSEMBLY);
204: MatAssemblyEnd(Amat,MAT_FINAL_ASSEMBLY);
205: }
206: }
207: } else {
208: if (ts->rhsjacobianp) {
209: TSComputeRHSJacobianP(ts,t,U,ts->Jacprhs);
210: }
211: if (ts->Jacprhs == Amat) { /* No IJacobian, so we only have the RHS matrix */
212: MatScale(Amat,-1);
213: } else if (ts->Jacprhs) { /* Both IJacobian and RHSJacobian */
214: MatStructure axpy = DIFFERENT_NONZERO_PATTERN;
215: if (!ts->ijacobianp) { /* No IJacobianp provided, but we have a separate RHS matrix */
216: MatZeroEntries(Amat);
217: }
218: MatAXPY(Amat,-1,ts->Jacprhs,axpy);
219: }
220: }
221: PetscLogEventEnd(TS_JacobianPEval,ts,U,Amat,0);
222: return(0);
223: }
225: /*@C
226: TSSetCostIntegrand - Sets the routine for evaluating the integral term in one or more cost functions
228: Logically Collective on TS
230: Input Parameters:
231: + ts - the TS context obtained from TSCreate()
232: . numcost - number of gradients to be computed, this is the number of cost functions
233: . costintegral - vector that stores the integral values
234: . rf - routine for evaluating the integrand function
235: . drduf - function that computes the gradients of the r's with respect to u
236: . 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)
237: . fwd - flag indicating whether to evaluate cost integral in the forward run or the adjoint run
238: - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
240: Calling sequence of rf:
241: $ PetscErrorCode rf(TS ts,PetscReal t,Vec U,Vec F,void *ctx);
243: Calling sequence of drduf:
244: $ PetscErroCode drduf(TS ts,PetscReal t,Vec U,Vec *dRdU,void *ctx);
246: Calling sequence of drdpf:
247: $ PetscErroCode drdpf(TS ts,PetscReal t,Vec U,Vec *dRdP,void *ctx);
249: Level: deprecated
251: Notes:
252: For optimization there is usually a single cost function (numcost = 1). For sensitivities there may be multiple cost functions
254: .seealso: TSSetRHSJacobianP(), TSGetCostGradients(), TSSetCostGradients()
255: @*/
256: PetscErrorCode TSSetCostIntegrand(TS ts,PetscInt numcost,Vec costintegral,PetscErrorCode (*rf)(TS,PetscReal,Vec,Vec,void*),
257: PetscErrorCode (*drduf)(TS,PetscReal,Vec,Vec*,void*),
258: PetscErrorCode (*drdpf)(TS,PetscReal,Vec,Vec*,void*),
259: PetscBool fwd,void *ctx)
260: {
266: 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()");
267: if (!ts->numcost) ts->numcost=numcost;
269: if (costintegral) {
270: PetscObjectReference((PetscObject)costintegral);
271: VecDestroy(&ts->vec_costintegral);
272: ts->vec_costintegral = costintegral;
273: } else {
274: if (!ts->vec_costintegral) { /* Create a seq vec if user does not provide one */
275: VecCreateSeq(PETSC_COMM_SELF,numcost,&ts->vec_costintegral);
276: } else {
277: VecSet(ts->vec_costintegral,0.0);
278: }
279: }
280: if (!ts->vec_costintegrand) {
281: VecDuplicate(ts->vec_costintegral,&ts->vec_costintegrand);
282: } else {
283: VecSet(ts->vec_costintegrand,0.0);
284: }
285: ts->costintegralfwd = fwd; /* Evaluate the cost integral in forward run if fwd is true */
286: ts->costintegrand = rf;
287: ts->costintegrandctx = ctx;
288: ts->drdufunction = drduf;
289: ts->drdpfunction = drdpf;
290: return(0);
291: }
293: /*@C
294: TSGetCostIntegral - Returns the values of the integral term in the cost functions.
295: It is valid to call the routine after a backward run.
297: Not Collective
299: Input Parameter:
300: . ts - the TS context obtained from TSCreate()
302: Output Parameter:
303: . v - the vector containing the integrals for each cost function
305: Level: intermediate
307: .seealso: TSSetCostIntegrand()
309: @*/
310: PetscErrorCode TSGetCostIntegral(TS ts,Vec *v)
311: {
312: TS quadts;
318: TSGetQuadratureTS(ts,NULL,&quadts);
319: *v = quadts->vec_sol;
320: return(0);
321: }
323: /*@C
324: TSComputeCostIntegrand - Evaluates the integral function in the cost functions.
326: Input Parameters:
327: + ts - the TS context
328: . t - current time
329: - U - state vector, i.e. current solution
331: Output Parameter:
332: . Q - vector of size numcost to hold the outputs
334: Notes:
335: Most users should not need to explicitly call this routine, as it
336: is used internally within the sensitivity analysis context.
338: Level: deprecated
340: .seealso: TSSetCostIntegrand()
341: @*/
342: PetscErrorCode TSComputeCostIntegrand(TS ts,PetscReal t,Vec U,Vec Q)
343: {
351: PetscLogEventBegin(TS_FunctionEval,ts,U,Q,0);
352: if (ts->costintegrand) {
353: PetscStackPush("TS user integrand in the cost function");
354: (*ts->costintegrand)(ts,t,U,Q,ts->costintegrandctx);
355: PetscStackPop;
356: } else {
357: VecZeroEntries(Q);
358: }
360: PetscLogEventEnd(TS_FunctionEval,ts,U,Q,0);
361: return(0);
362: }
364: /*@C
365: TSComputeDRDUFunction - Deprecated, use TSGetQuadratureTS() then TSComputeRHSJacobian()
367: Level: deprecated
369: @*/
370: PetscErrorCode TSComputeDRDUFunction(TS ts,PetscReal t,Vec U,Vec *DRDU)
371: {
375: if (!DRDU) return(0);
379: PetscStackPush("TS user DRDU function for sensitivity analysis");
380: (*ts->drdufunction)(ts,t,U,DRDU,ts->costintegrandctx);
381: PetscStackPop;
382: return(0);
383: }
385: /*@C
386: TSComputeDRDPFunction - Deprecated, use TSGetQuadratureTS() then TSComputeRHSJacobianP()
388: Level: deprecated
390: @*/
391: PetscErrorCode TSComputeDRDPFunction(TS ts,PetscReal t,Vec U,Vec *DRDP)
392: {
396: if (!DRDP) return(0);
400: PetscStackPush("TS user DRDP function for sensitivity analysis");
401: (*ts->drdpfunction)(ts,t,U,DRDP,ts->costintegrandctx);
402: PetscStackPop;
403: return(0);
404: }
406: /*@C
407: 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.
409: Logically Collective on TS
411: Input Parameters:
412: + ts - TS context obtained from TSCreate()
413: . ihp1 - an array of vectors storing the result of vector-Hessian-vector product for F_UU
414: . hessianproductfunc1 - vector-Hessian-vector product function for F_UU
415: . ihp2 - an array of vectors storing the result of vector-Hessian-vector product for F_UP
416: . hessianproductfunc2 - vector-Hessian-vector product function for F_UP
417: . ihp3 - an array of vectors storing the result of vector-Hessian-vector product for F_PU
418: . hessianproductfunc3 - vector-Hessian-vector product function for F_PU
419: . ihp4 - an array of vectors storing the result of vector-Hessian-vector product for F_PP
420: . hessianproductfunc4 - vector-Hessian-vector product function for F_PP
422: Calling sequence of ihessianproductfunc:
423: $ ihessianproductfunc (TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx);
424: + t - current timestep
425: . U - input vector (current ODE solution)
426: . Vl - an array of input vectors to be left-multiplied with the Hessian
427: . Vr - input vector to be right-multiplied with the Hessian
428: . VHV - an array of output vectors for vector-Hessian-vector product
429: - ctx - [optional] user-defined function context
431: Level: intermediate
433: Notes:
434: The first Hessian function and the working array are required.
435: As an example to implement the callback functions, the second callback function calculates the vector-Hessian-vector product
436: $ Vl_n^T*F_UP*Vr
437: where the vector Vl_n (n-th element in the array Vl) and Vr are of size N and M respectively, and the Hessian F_UP is of size N x N x M.
438: Each entry of F_UP corresponds to the derivative
439: $ F_UP[i][j][k] = \frac{\partial^2 F[i]}{\partial U[j] \partial P[k]}.
440: The result of the vector-Hessian-vector product for Vl_n needs to be stored in vector VHV_n with the j-th entry being
441: $ VHV_n[j] = \sum_i \sum_k {Vl_n[i] * F_UP[i][j][k] * Vr[k]}
442: If the cost function is a scalar, there will be only one vector in Vl and VHV.
444: .seealso:
445: @*/
446: PetscErrorCode TSSetIHessianProduct(TS ts,Vec *ihp1,PetscErrorCode (*ihessianproductfunc1)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*),
447: Vec *ihp2,PetscErrorCode (*ihessianproductfunc2)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*),
448: Vec *ihp3,PetscErrorCode (*ihessianproductfunc3)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*),
449: Vec *ihp4,PetscErrorCode (*ihessianproductfunc4)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*),
450: void *ctx)
451: {
456: ts->ihessianproductctx = ctx;
457: if (ihp1) ts->vecs_fuu = ihp1;
458: if (ihp2) ts->vecs_fup = ihp2;
459: if (ihp3) ts->vecs_fpu = ihp3;
460: if (ihp4) ts->vecs_fpp = ihp4;
461: ts->ihessianproduct_fuu = ihessianproductfunc1;
462: ts->ihessianproduct_fup = ihessianproductfunc2;
463: ts->ihessianproduct_fpu = ihessianproductfunc3;
464: ts->ihessianproduct_fpp = ihessianproductfunc4;
465: return(0);
466: }
468: /*@C
469: TSComputeIHessianProductFunctionUU - Runs the user-defined vector-Hessian-vector product function for Fuu.
471: Collective on TS
473: Input Parameters:
474: . ts - The TS context obtained from TSCreate()
476: Notes:
477: TSComputeIHessianProductFunctionUU() is typically used for sensitivity implementation,
478: so most users would not generally call this routine themselves.
480: Level: developer
482: .seealso: TSSetIHessianProduct()
483: @*/
484: PetscErrorCode TSComputeIHessianProductFunctionUU(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV)
485: {
489: if (!VHV) return(0);
493: if (ts->ihessianproduct_fuu) {
494: PetscStackPush("TS user IHessianProduct function 1 for sensitivity analysis");
495: (*ts->ihessianproduct_fuu)(ts,t,U,Vl,Vr,VHV,ts->ihessianproductctx);
496: PetscStackPop;
497: }
498: /* does not consider IMEX for now, so either IHessian or RHSHessian will be calculated, using the same output VHV */
499: if (ts->rhshessianproduct_guu) {
500: PetscInt nadj;
501: TSComputeRHSHessianProductFunctionUU(ts,t,U,Vl,Vr,VHV);
502: for (nadj=0; nadj<ts->numcost; nadj++) {
503: VecScale(VHV[nadj],-1);
504: }
505: }
506: return(0);
507: }
509: /*@C
510: TSComputeIHessianProductFunctionUP - Runs the user-defined vector-Hessian-vector product function for Fup.
512: Collective on TS
514: Input Parameters:
515: . ts - The TS context obtained from TSCreate()
517: Notes:
518: TSComputeIHessianProductFunctionUP() is typically used for sensitivity implementation,
519: so most users would not generally call this routine themselves.
521: Level: developer
523: .seealso: TSSetIHessianProduct()
524: @*/
525: PetscErrorCode TSComputeIHessianProductFunctionUP(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV)
526: {
530: if (!VHV) return(0);
534: if (ts->ihessianproduct_fup) {
535: PetscStackPush("TS user IHessianProduct function 2 for sensitivity analysis");
536: (*ts->ihessianproduct_fup)(ts,t,U,Vl,Vr,VHV,ts->ihessianproductctx);
537: PetscStackPop;
538: }
539: /* does not consider IMEX for now, so either IHessian or RHSHessian will be calculated, using the same output VHV */
540: if (ts->rhshessianproduct_gup) {
541: PetscInt nadj;
542: TSComputeRHSHessianProductFunctionUP(ts,t,U,Vl,Vr,VHV);
543: for (nadj=0; nadj<ts->numcost; nadj++) {
544: VecScale(VHV[nadj],-1);
545: }
546: }
547: return(0);
548: }
550: /*@C
551: TSComputeIHessianProductFunctionPU - Runs the user-defined vector-Hessian-vector product function for Fpu.
553: Collective on TS
555: Input Parameters:
556: . ts - The TS context obtained from TSCreate()
558: Notes:
559: TSComputeIHessianProductFunctionPU() is typically used for sensitivity implementation,
560: so most users would not generally call this routine themselves.
562: Level: developer
564: .seealso: TSSetIHessianProduct()
565: @*/
566: PetscErrorCode TSComputeIHessianProductFunctionPU(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV)
567: {
571: if (!VHV) return(0);
575: if (ts->ihessianproduct_fpu) {
576: PetscStackPush("TS user IHessianProduct function 3 for sensitivity analysis");
577: (*ts->ihessianproduct_fpu)(ts,t,U,Vl,Vr,VHV,ts->ihessianproductctx);
578: PetscStackPop;
579: }
580: /* does not consider IMEX for now, so either IHessian or RHSHessian will be calculated, using the same output VHV */
581: if (ts->rhshessianproduct_gpu) {
582: PetscInt nadj;
583: TSComputeRHSHessianProductFunctionPU(ts,t,U,Vl,Vr,VHV);
584: for (nadj=0; nadj<ts->numcost; nadj++) {
585: VecScale(VHV[nadj],-1);
586: }
587: }
588: return(0);
589: }
591: /*@C
592: TSComputeIHessianProductFunctionPP - Runs the user-defined vector-Hessian-vector product function for Fpp.
594: Collective on TS
596: Input Parameters:
597: . ts - The TS context obtained from TSCreate()
599: Notes:
600: TSComputeIHessianProductFunctionPP() is typically used for sensitivity implementation,
601: so most users would not generally call this routine themselves.
603: Level: developer
605: .seealso: TSSetIHessianProduct()
606: @*/
607: PetscErrorCode TSComputeIHessianProductFunctionPP(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV)
608: {
612: if (!VHV) return(0);
616: if (ts->ihessianproduct_fpp) {
617: PetscStackPush("TS user IHessianProduct function 3 for sensitivity analysis");
618: (*ts->ihessianproduct_fpp)(ts,t,U,Vl,Vr,VHV,ts->ihessianproductctx);
619: PetscStackPop;
620: }
621: /* does not consider IMEX for now, so either IHessian or RHSHessian will be calculated, using the same output VHV */
622: if (ts->rhshessianproduct_gpp) {
623: PetscInt nadj;
624: TSComputeRHSHessianProductFunctionPP(ts,t,U,Vl,Vr,VHV);
625: for (nadj=0; nadj<ts->numcost; nadj++) {
626: VecScale(VHV[nadj],-1);
627: }
628: }
629: return(0);
630: }
632: /*@C
633: TSSetRHSHessianProduct - Sets the function that computes the vector-Hessian-vector product. The Hessian is the second-order derivative of G (RHSFunction) w.r.t. the state variable.
635: Logically Collective on TS
637: Input Parameters:
638: + ts - TS context obtained from TSCreate()
639: . rhshp1 - an array of vectors storing the result of vector-Hessian-vector product for G_UU
640: . hessianproductfunc1 - vector-Hessian-vector product function for G_UU
641: . rhshp2 - an array of vectors storing the result of vector-Hessian-vector product for G_UP
642: . hessianproductfunc2 - vector-Hessian-vector product function for G_UP
643: . rhshp3 - an array of vectors storing the result of vector-Hessian-vector product for G_PU
644: . hessianproductfunc3 - vector-Hessian-vector product function for G_PU
645: . rhshp4 - an array of vectors storing the result of vector-Hessian-vector product for G_PP
646: . hessianproductfunc4 - vector-Hessian-vector product function for G_PP
648: Calling sequence of ihessianproductfunc:
649: $ rhshessianproductfunc (TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx);
650: + t - current timestep
651: . U - input vector (current ODE solution)
652: . Vl - an array of input vectors to be left-multiplied with the Hessian
653: . Vr - input vector to be right-multiplied with the Hessian
654: . VHV - an array of output vectors for vector-Hessian-vector product
655: - ctx - [optional] user-defined function context
657: Level: intermediate
659: Notes:
660: The first Hessian function and the working array are required.
661: As an example to implement the callback functions, the second callback function calculates the vector-Hessian-vector product
662: $ Vl_n^T*G_UP*Vr
663: where the vector Vl_n (n-th element in the array Vl) and Vr are of size N and M respectively, and the Hessian G_UP is of size N x N x M.
664: Each entry of G_UP corresponds to the derivative
665: $ G_UP[i][j][k] = \frac{\partial^2 G[i]}{\partial U[j] \partial P[k]}.
666: The result of the vector-Hessian-vector product for Vl_n needs to be stored in vector VHV_n with j-th entry being
667: $ VHV_n[j] = \sum_i \sum_k {Vl_n[i] * G_UP[i][j][k] * Vr[k]}
668: If the cost function is a scalar, there will be only one vector in Vl and VHV.
670: .seealso:
671: @*/
672: PetscErrorCode TSSetRHSHessianProduct(TS ts,Vec *rhshp1,PetscErrorCode (*rhshessianproductfunc1)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*),
673: Vec *rhshp2,PetscErrorCode (*rhshessianproductfunc2)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*),
674: Vec *rhshp3,PetscErrorCode (*rhshessianproductfunc3)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*),
675: Vec *rhshp4,PetscErrorCode (*rhshessianproductfunc4)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*),
676: void *ctx)
677: {
682: ts->rhshessianproductctx = ctx;
683: if (rhshp1) ts->vecs_guu = rhshp1;
684: if (rhshp2) ts->vecs_gup = rhshp2;
685: if (rhshp3) ts->vecs_gpu = rhshp3;
686: if (rhshp4) ts->vecs_gpp = rhshp4;
687: ts->rhshessianproduct_guu = rhshessianproductfunc1;
688: ts->rhshessianproduct_gup = rhshessianproductfunc2;
689: ts->rhshessianproduct_gpu = rhshessianproductfunc3;
690: ts->rhshessianproduct_gpp = rhshessianproductfunc4;
691: return(0);
692: }
694: /*@C
695: TSComputeRHSHessianProductFunctionUU - Runs the user-defined vector-Hessian-vector product function for Guu.
697: Collective on TS
699: Input Parameters:
700: . ts - The TS context obtained from TSCreate()
702: Notes:
703: TSComputeRHSHessianProductFunctionUU() is typically used for sensitivity implementation,
704: so most users would not generally call this routine themselves.
706: Level: developer
708: .seealso: TSSetRHSHessianProduct()
709: @*/
710: PetscErrorCode TSComputeRHSHessianProductFunctionUU(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV)
711: {
715: if (!VHV) return(0);
719: PetscStackPush("TS user RHSHessianProduct function 1 for sensitivity analysis");
720: (*ts->rhshessianproduct_guu)(ts,t,U,Vl,Vr,VHV,ts->rhshessianproductctx);
721: PetscStackPop;
722: return(0);
723: }
725: /*@C
726: TSComputeRHSHessianProductFunctionUP - Runs the user-defined vector-Hessian-vector product function for Gup.
728: Collective on TS
730: Input Parameters:
731: . ts - The TS context obtained from TSCreate()
733: Notes:
734: TSComputeRHSHessianProductFunctionUP() is typically used for sensitivity implementation,
735: so most users would not generally call this routine themselves.
737: Level: developer
739: .seealso: TSSetRHSHessianProduct()
740: @*/
741: PetscErrorCode TSComputeRHSHessianProductFunctionUP(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV)
742: {
746: if (!VHV) return(0);
750: PetscStackPush("TS user RHSHessianProduct function 2 for sensitivity analysis");
751: (*ts->rhshessianproduct_gup)(ts,t,U,Vl,Vr,VHV,ts->rhshessianproductctx);
752: PetscStackPop;
753: return(0);
754: }
756: /*@C
757: TSComputeRHSHessianProductFunctionPU - Runs the user-defined vector-Hessian-vector product function for Gpu.
759: Collective on TS
761: Input Parameters:
762: . ts - The TS context obtained from TSCreate()
764: Notes:
765: TSComputeRHSHessianProductFunctionPU() is typically used for sensitivity implementation,
766: so most users would not generally call this routine themselves.
768: Level: developer
770: .seealso: TSSetRHSHessianProduct()
771: @*/
772: PetscErrorCode TSComputeRHSHessianProductFunctionPU(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV)
773: {
777: if (!VHV) return(0);
781: PetscStackPush("TS user RHSHessianProduct function 3 for sensitivity analysis");
782: (*ts->rhshessianproduct_gpu)(ts,t,U,Vl,Vr,VHV,ts->rhshessianproductctx);
783: PetscStackPop;
784: return(0);
785: }
787: /*@C
788: TSComputeRHSHessianProductFunctionPP - Runs the user-defined vector-Hessian-vector product function for Gpp.
790: Collective on TS
792: Input Parameters:
793: . ts - The TS context obtained from TSCreate()
795: Notes:
796: TSComputeRHSHessianProductFunctionPP() is typically used for sensitivity implementation,
797: so most users would not generally call this routine themselves.
799: Level: developer
801: .seealso: TSSetRHSHessianProduct()
802: @*/
803: PetscErrorCode TSComputeRHSHessianProductFunctionPP(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV)
804: {
808: if (!VHV) return(0);
812: PetscStackPush("TS user RHSHessianProduct function 3 for sensitivity analysis");
813: (*ts->rhshessianproduct_gpp)(ts,t,U,Vl,Vr,VHV,ts->rhshessianproductctx);
814: PetscStackPop;
815: return(0);
816: }
818: /* --------------------------- Adjoint sensitivity ---------------------------*/
820: /*@
821: TSSetCostGradients - Sets the initial value of the gradients of the cost function w.r.t. initial values and w.r.t. the problem parameters
822: for use by the TSAdjoint routines.
824: Logically Collective on TS
826: Input Parameters:
827: + ts - the TS context obtained from TSCreate()
828: . numcost - number of gradients to be computed, this is the number of cost functions
829: . 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
830: - mu - gradients with respect to the parameters, the number of entries in these vectors is the same as the number of parameters
832: Level: beginner
834: Notes:
835: the entries in these vectors must be correctly initialized with the values lamda_i = df/dy|finaltime mu_i = df/dp|finaltime
837: After TSAdjointSolve() is called the lamba and the mu contain the computed sensitivities
839: .seealso TSGetCostGradients()
840: @*/
841: PetscErrorCode TSSetCostGradients(TS ts,PetscInt numcost,Vec *lambda,Vec *mu)
842: {
846: ts->vecs_sensi = lambda;
847: ts->vecs_sensip = mu;
848: 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");
849: ts->numcost = numcost;
850: return(0);
851: }
853: /*@
854: TSGetCostGradients - Returns the gradients from the TSAdjointSolve()
856: Not Collective, but Vec returned is parallel if TS is parallel
858: Input Parameter:
859: . ts - the TS context obtained from TSCreate()
861: Output Parameter:
862: + lambda - vectors containing the gradients of the cost functions with respect to the ODE/DAE solution variables
863: - mu - vectors containing the gradients of the cost functions with respect to the problem parameters
865: Level: intermediate
867: .seealso: TSSetCostGradients()
868: @*/
869: PetscErrorCode TSGetCostGradients(TS ts,PetscInt *numcost,Vec **lambda,Vec **mu)
870: {
873: if (numcost) *numcost = ts->numcost;
874: if (lambda) *lambda = ts->vecs_sensi;
875: if (mu) *mu = ts->vecs_sensip;
876: return(0);
877: }
879: /*@
880: 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
881: for use by the TSAdjoint routines.
883: Logically Collective on TS
885: Input Parameters:
886: + ts - the TS context obtained from TSCreate()
887: . numcost - number of cost functions
888: . 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
889: . mu2 - Hessian-vector product with respect to the parameters, the number of entries in these vectors is the same as the number of parameters
890: - dir - the direction vector that are multiplied with the Hessian of the cost functions
892: Level: beginner
894: Notes: Hessian of the cost function is completely different from Hessian of the ODE/DAE system
896: For second-order adjoint, one needs to call this function and then TSAdjointSetForward() before TSSolve().
898: 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.
900: Passing NULL for lambda2 disables the second-order calculation.
901: .seealso: TSAdjointSetForward()
902: @*/
903: PetscErrorCode TSSetCostHessianProducts(TS ts,PetscInt numcost,Vec *lambda2,Vec *mu2,Vec dir)
904: {
907: 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");
908: ts->numcost = numcost;
909: ts->vecs_sensi2 = lambda2;
910: ts->vecs_sensi2p = mu2;
911: ts->vec_dir = dir;
912: return(0);
913: }
915: /*@
916: TSGetCostHessianProducts - Returns the gradients from the TSAdjointSolve()
918: Not Collective, but Vec returned is parallel if TS is parallel
920: Input Parameter:
921: . ts - the TS context obtained from TSCreate()
923: Output Parameter:
924: + numcost - number of cost functions
925: . 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
926: . mu2 - Hessian-vector product with respect to the parameters, the number of entries in these vectors is the same as the number of parameters
927: - dir - the direction vector that are multiplied with the Hessian of the cost functions
929: Level: intermediate
931: .seealso: TSSetCostHessianProducts()
932: @*/
933: PetscErrorCode TSGetCostHessianProducts(TS ts,PetscInt *numcost,Vec **lambda2,Vec **mu2, Vec *dir)
934: {
937: if (numcost) *numcost = ts->numcost;
938: if (lambda2) *lambda2 = ts->vecs_sensi2;
939: if (mu2) *mu2 = ts->vecs_sensi2p;
940: if (dir) *dir = ts->vec_dir;
941: return(0);
942: }
944: /*@
945: TSAdjointSetForward - Trigger the tangent linear solver and initialize the forward sensitivities
947: Logically Collective on TS
949: Input Parameters:
950: + ts - the TS context obtained from TSCreate()
951: - didp - the derivative of initial values w.r.t. parameters
953: Level: intermediate
955: Notes: When computing sensitivies w.r.t. initial condition, set didp to NULL so that the solver will take it as an identity matrix mathematically. TSAdjoint does not reset the tangent linear solver automatically, TSAdjointResetForward() should be called to reset the tangent linear solver.
957: .seealso: TSSetCostHessianProducts(), TSAdjointResetForward()
958: @*/
959: PetscErrorCode TSAdjointSetForward(TS ts,Mat didp)
960: {
961: Mat A;
962: Vec sp;
963: PetscScalar *xarr;
964: PetscInt lsize;
968: ts->forward_solve = PETSC_TRUE; /* turn on tangent linear mode */
969: if (!ts->vecs_sensi2) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must call TSSetCostHessianProducts() first");
970: if (!ts->vec_dir) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Directional vector is missing. Call TSSetCostHessianProducts() to set it.");
971: /* create a single-column dense matrix */
972: VecGetLocalSize(ts->vec_sol,&lsize);
973: MatCreateDense(PetscObjectComm((PetscObject)ts),lsize,PETSC_DECIDE,PETSC_DECIDE,1,NULL,&A);
975: VecDuplicate(ts->vec_sol,&sp);
976: MatDenseGetColumn(A,0,&xarr);
977: VecPlaceArray(sp,xarr);
978: if (ts->vecs_sensi2p) { /* tangent linear variable initialized as 2*dIdP*dir */
979: if (didp) {
980: MatMult(didp,ts->vec_dir,sp);
981: VecScale(sp,2.);
982: } else {
983: VecZeroEntries(sp);
984: }
985: } else { /* tangent linear variable initialized as dir */
986: VecCopy(ts->vec_dir,sp);
987: }
988: VecResetArray(sp);
989: MatDenseRestoreColumn(A,&xarr);
990: VecDestroy(&sp);
992: TSForwardSetInitialSensitivities(ts,A); /* if didp is NULL, identity matrix is assumed */
994: MatDestroy(&A);
995: return(0);
996: }
998: /*@
999: TSAdjointResetForward - Reset the tangent linear solver and destroy the tangent linear context
1001: Logically Collective on TS
1003: Input Parameters:
1004: . ts - the TS context obtained from TSCreate()
1006: Level: intermediate
1008: .seealso: TSAdjointSetForward()
1009: @*/
1010: PetscErrorCode TSAdjointResetForward(TS ts)
1011: {
1015: ts->forward_solve = PETSC_FALSE; /* turn off tangent linear mode */
1016: TSForwardReset(ts);
1017: return(0);
1018: }
1020: /*@
1021: TSAdjointSetUp - Sets up the internal data structures for the later use
1022: of an adjoint solver
1024: Collective on TS
1026: Input Parameter:
1027: . ts - the TS context obtained from TSCreate()
1029: Level: advanced
1031: .seealso: TSCreate(), TSAdjointStep(), TSSetCostGradients()
1032: @*/
1033: PetscErrorCode TSAdjointSetUp(TS ts)
1034: {
1035: TSTrajectory tj;
1036: PetscBool match;
1037: PetscErrorCode ierr;
1041: if (ts->adjointsetupcalled) return(0);
1042: if (!ts->vecs_sensi) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetCostGradients() first");
1043: if (ts->vecs_sensip && !ts->Jacp && !ts->Jacprhs) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetRHSJacobianP() or TSSetIJacobianP() first");
1044: TSGetTrajectory(ts,&tj);
1045: PetscObjectTypeCompare((PetscObject)tj,TSTRAJECTORYBASIC,&match);
1046: if (match) {
1047: PetscBool solution_only;
1048: TSTrajectoryGetSolutionOnly(tj,&solution_only);
1049: if (solution_only) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"TSAdjoint cannot use the solution-only mode when choosing the Basic TSTrajectory type. Turn it off with -ts_trajectory_solution_only 0");
1050: }
1051: TSTrajectorySetUseHistory(tj,PETSC_FALSE); /* not use TSHistory */
1053: if (ts->quadraturets) { /* if there is integral in the cost function */
1054: VecDuplicate(ts->vecs_sensi[0],&ts->vec_drdu_col);
1055: if (ts->vecs_sensip){
1056: VecDuplicate(ts->vecs_sensip[0],&ts->vec_drdp_col);
1057: }
1058: }
1060: if (ts->ops->adjointsetup) {
1061: (*ts->ops->adjointsetup)(ts);
1062: }
1063: ts->adjointsetupcalled = PETSC_TRUE;
1064: return(0);
1065: }
1067: /*@
1068: TSAdjointReset - Resets a TSAdjoint context and removes any allocated Vecs and Mats.
1070: Collective on TS
1072: Input Parameter:
1073: . ts - the TS context obtained from TSCreate()
1075: Level: beginner
1077: .seealso: TSCreate(), TSAdjointSetUp(), TSADestroy()
1078: @*/
1079: PetscErrorCode TSAdjointReset(TS ts)
1080: {
1085: if (ts->ops->adjointreset) {
1086: (*ts->ops->adjointreset)(ts);
1087: }
1088: if (ts->quadraturets) { /* if there is integral in the cost function */
1089: VecDestroy(&ts->vec_drdu_col);
1090: if (ts->vecs_sensip){
1091: VecDestroy(&ts->vec_drdp_col);
1092: }
1093: }
1094: ts->vecs_sensi = NULL;
1095: ts->vecs_sensip = NULL;
1096: ts->vecs_sensi2 = NULL;
1097: ts->vecs_sensi2p = NULL;
1098: ts->vec_dir = NULL;
1099: ts->adjointsetupcalled = PETSC_FALSE;
1100: return(0);
1101: }
1103: /*@
1104: TSAdjointSetSteps - Sets the number of steps the adjoint solver should take backward in time
1106: Logically Collective on TS
1108: Input Parameters:
1109: + ts - the TS context obtained from TSCreate()
1110: - steps - number of steps to use
1112: Level: intermediate
1114: Notes:
1115: Normally one does not call this and TSAdjointSolve() integrates back to the original timestep. One can call this
1116: so as to integrate back to less than the original timestep
1118: .seealso: TSSetExactFinalTime()
1119: @*/
1120: PetscErrorCode TSAdjointSetSteps(TS ts,PetscInt steps)
1121: {
1125: if (steps < 0) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Cannot step back a negative number of steps");
1126: if (steps > ts->steps) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Cannot step back more than the total number of forward steps");
1127: ts->adjoint_max_steps = steps;
1128: return(0);
1129: }
1131: /*@C
1132: TSAdjointSetRHSJacobian - Deprecated, use TSSetRHSJacobianP()
1134: Level: deprecated
1136: @*/
1137: PetscErrorCode TSAdjointSetRHSJacobian(TS ts,Mat Amat,PetscErrorCode (*func)(TS,PetscReal,Vec,Mat,void*),void *ctx)
1138: {
1145: ts->rhsjacobianp = func;
1146: ts->rhsjacobianpctx = ctx;
1147: if(Amat) {
1148: PetscObjectReference((PetscObject)Amat);
1149: MatDestroy(&ts->Jacp);
1150: ts->Jacp = Amat;
1151: }
1152: return(0);
1153: }
1155: /*@C
1156: TSAdjointComputeRHSJacobian - Deprecated, use TSComputeRHSJacobianP()
1158: Level: deprecated
1160: @*/
1161: PetscErrorCode TSAdjointComputeRHSJacobian(TS ts,PetscReal t,Vec U,Mat Amat)
1162: {
1170: PetscStackPush("TS user JacobianP function for sensitivity analysis");
1171: (*ts->rhsjacobianp)(ts,t,U,Amat,ts->rhsjacobianpctx);
1172: PetscStackPop;
1173: return(0);
1174: }
1176: /*@
1177: TSAdjointComputeDRDYFunction - Deprecated, use TSGetQuadratureTS() then TSComputeRHSJacobian()
1179: Level: deprecated
1181: @*/
1182: PetscErrorCode TSAdjointComputeDRDYFunction(TS ts,PetscReal t,Vec U,Vec *DRDU)
1183: {
1190: PetscStackPush("TS user DRDY function for sensitivity analysis");
1191: (*ts->drdufunction)(ts,t,U,DRDU,ts->costintegrandctx);
1192: PetscStackPop;
1193: return(0);
1194: }
1196: /*@
1197: TSAdjointComputeDRDPFunction - Deprecated, use TSGetQuadratureTS() then TSComputeRHSJacobianP()
1199: Level: deprecated
1201: @*/
1202: PetscErrorCode TSAdjointComputeDRDPFunction(TS ts,PetscReal t,Vec U,Vec *DRDP)
1203: {
1210: PetscStackPush("TS user DRDP function for sensitivity analysis");
1211: (*ts->drdpfunction)(ts,t,U,DRDP,ts->costintegrandctx);
1212: PetscStackPop;
1213: return(0);
1214: }
1216: /*@C
1217: TSAdjointMonitorSensi - monitors the first lambda sensitivity
1219: Level: intermediate
1221: .seealso: TSAdjointMonitorSet()
1222: @*/
1223: PetscErrorCode TSAdjointMonitorSensi(TS ts,PetscInt step,PetscReal ptime,Vec v,PetscInt numcost,Vec *lambda,Vec *mu,PetscViewerAndFormat *vf)
1224: {
1226: PetscViewer viewer = vf->viewer;
1230: PetscViewerPushFormat(viewer,vf->format);
1231: VecView(lambda[0],viewer);
1232: PetscViewerPopFormat(viewer);
1233: return(0);
1234: }
1236: /*@C
1237: TSAdjointMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user
1239: Collective on TS
1241: Input Parameters:
1242: + ts - TS object you wish to monitor
1243: . name - the monitor type one is seeking
1244: . help - message indicating what monitoring is done
1245: . manual - manual page for the monitor
1246: . monitor - the monitor function
1247: - 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
1249: Level: developer
1251: .seealso: PetscOptionsGetViewer(), PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(),
1252: PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
1253: PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
1254: PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
1255: PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1256: PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1257: PetscOptionsFList(), PetscOptionsEList()
1258: @*/
1259: 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*))
1260: {
1261: PetscErrorCode ierr;
1262: PetscViewer viewer;
1263: PetscViewerFormat format;
1264: PetscBool flg;
1267: PetscOptionsGetViewer(PetscObjectComm((PetscObject)ts),((PetscObject) ts)->options,((PetscObject)ts)->prefix,name,&viewer,&format,&flg);
1268: if (flg) {
1269: PetscViewerAndFormat *vf;
1270: PetscViewerAndFormatCreate(viewer,format,&vf);
1271: PetscObjectDereference((PetscObject)viewer);
1272: if (monitorsetup) {
1273: (*monitorsetup)(ts,vf);
1274: }
1275: TSAdjointMonitorSet(ts,(PetscErrorCode (*)(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*))monitor,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy);
1276: }
1277: return(0);
1278: }
1280: /*@C
1281: TSAdjointMonitorSet - Sets an ADDITIONAL function that is to be used at every
1282: timestep to display the iteration's progress.
1284: Logically Collective on TS
1286: Input Parameters:
1287: + ts - the TS context obtained from TSCreate()
1288: . adjointmonitor - monitoring routine
1289: . adjointmctx - [optional] user-defined context for private data for the
1290: monitor routine (use NULL if no context is desired)
1291: - adjointmonitordestroy - [optional] routine that frees monitor context
1292: (may be NULL)
1294: Calling sequence of monitor:
1295: $ int adjointmonitor(TS ts,PetscInt steps,PetscReal time,Vec u,PetscInt numcost,Vec *lambda, Vec *mu,void *adjointmctx)
1297: + ts - the TS context
1298: . 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
1299: been interpolated to)
1300: . time - current time
1301: . u - current iterate
1302: . numcost - number of cost functionos
1303: . lambda - sensitivities to initial conditions
1304: . mu - sensitivities to parameters
1305: - adjointmctx - [optional] adjoint monitoring context
1307: Notes:
1308: This routine adds an additional monitor to the list of monitors that
1309: already has been loaded.
1311: Fortran Notes:
1312: Only a single monitor function can be set for each TS object
1314: Level: intermediate
1316: .seealso: TSAdjointMonitorCancel()
1317: @*/
1318: PetscErrorCode TSAdjointMonitorSet(TS ts,PetscErrorCode (*adjointmonitor)(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*),void *adjointmctx,PetscErrorCode (*adjointmdestroy)(void**))
1319: {
1321: PetscInt i;
1322: PetscBool identical;
1326: for (i=0; i<ts->numbermonitors;i++) {
1327: PetscMonitorCompare((PetscErrorCode (*)(void))adjointmonitor,adjointmctx,adjointmdestroy,(PetscErrorCode (*)(void))ts->adjointmonitor[i],ts->adjointmonitorcontext[i],ts->adjointmonitordestroy[i],&identical);
1328: if (identical) return(0);
1329: }
1330: if (ts->numberadjointmonitors >= MAXTSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many adjoint monitors set");
1331: ts->adjointmonitor[ts->numberadjointmonitors] = adjointmonitor;
1332: ts->adjointmonitordestroy[ts->numberadjointmonitors] = adjointmdestroy;
1333: ts->adjointmonitorcontext[ts->numberadjointmonitors++] = (void*)adjointmctx;
1334: return(0);
1335: }
1337: /*@C
1338: TSAdjointMonitorCancel - Clears all the adjoint monitors that have been set on a time-step object.
1340: Logically Collective on TS
1342: Input Parameters:
1343: . ts - the TS context obtained from TSCreate()
1345: Notes:
1346: There is no way to remove a single, specific monitor.
1348: Level: intermediate
1350: .seealso: TSAdjointMonitorSet()
1351: @*/
1352: PetscErrorCode TSAdjointMonitorCancel(TS ts)
1353: {
1355: PetscInt i;
1359: for (i=0; i<ts->numberadjointmonitors; i++) {
1360: if (ts->adjointmonitordestroy[i]) {
1361: (*ts->adjointmonitordestroy[i])(&ts->adjointmonitorcontext[i]);
1362: }
1363: }
1364: ts->numberadjointmonitors = 0;
1365: return(0);
1366: }
1368: /*@C
1369: TSAdjointMonitorDefault - the default monitor of adjoint computations
1371: Level: intermediate
1373: .seealso: TSAdjointMonitorSet()
1374: @*/
1375: PetscErrorCode TSAdjointMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,PetscInt numcost,Vec *lambda,Vec *mu,PetscViewerAndFormat *vf)
1376: {
1378: PetscViewer viewer = vf->viewer;
1382: PetscViewerPushFormat(viewer,vf->format);
1383: PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);
1384: PetscViewerASCIIPrintf(viewer,"%D TS dt %g time %g%s",step,(double)ts->time_step,(double)ptime,ts->steprollback ? " (r)\n" : "\n");
1385: PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);
1386: PetscViewerPopFormat(viewer);
1387: return(0);
1388: }
1390: /*@C
1391: TSAdjointMonitorDrawSensi - Monitors progress of the adjoint TS solvers by calling
1392: VecView() for the sensitivities to initial states at each timestep
1394: Collective on TS
1396: Input Parameters:
1397: + ts - the TS context
1398: . step - current time-step
1399: . ptime - current time
1400: . u - current state
1401: . numcost - number of cost functions
1402: . lambda - sensitivities to initial conditions
1403: . mu - sensitivities to parameters
1404: - dummy - either a viewer or NULL
1406: Level: intermediate
1408: .seealso: TSAdjointMonitorSet(), TSAdjointMonitorDefault(), VecView()
1409: @*/
1410: PetscErrorCode TSAdjointMonitorDrawSensi(TS ts,PetscInt step,PetscReal ptime,Vec u,PetscInt numcost,Vec *lambda,Vec *mu,void *dummy)
1411: {
1412: PetscErrorCode ierr;
1413: TSMonitorDrawCtx ictx = (TSMonitorDrawCtx)dummy;
1414: PetscDraw draw;
1415: PetscReal xl,yl,xr,yr,h;
1416: char time[32];
1419: if (!(((ictx->howoften > 0) && (!(step % ictx->howoften))) || ((ictx->howoften == -1) && ts->reason))) return(0);
1421: VecView(lambda[0],ictx->viewer);
1422: PetscViewerDrawGetDraw(ictx->viewer,0,&draw);
1423: PetscSNPrintf(time,32,"Timestep %d Time %g",(int)step,(double)ptime);
1424: PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);
1425: h = yl + .95*(yr - yl);
1426: PetscDrawStringCentered(draw,.5*(xl+xr),h,PETSC_DRAW_BLACK,time);
1427: PetscDrawFlush(draw);
1428: return(0);
1429: }
1431: /*
1432: TSAdjointSetFromOptions - Sets various TSAdjoint parameters from user options.
1434: Collective on TSAdjoint
1436: Input Parameter:
1437: . ts - the TS context
1439: Options Database Keys:
1440: + -ts_adjoint_solve <yes,no> After solving the ODE/DAE solve the adjoint problem (requires -ts_save_trajectory)
1441: . -ts_adjoint_monitor - print information at each adjoint time step
1442: - -ts_adjoint_monitor_draw_sensi - monitor the sensitivity of the first cost function wrt initial conditions (lambda[0]) graphically
1444: Level: developer
1446: Notes:
1447: This is not normally called directly by users
1449: .seealso: TSSetSaveTrajectory(), TSTrajectorySetUp()
1450: */
1451: PetscErrorCode TSAdjointSetFromOptions(PetscOptionItems *PetscOptionsObject,TS ts)
1452: {
1453: PetscBool tflg,opt;
1458: PetscOptionsHead(PetscOptionsObject,"TS Adjoint options");
1459: tflg = ts->adjoint_solve ? PETSC_TRUE : PETSC_FALSE;
1460: PetscOptionsBool("-ts_adjoint_solve","Solve the adjoint problem immediately after solving the forward problem","",tflg,&tflg,&opt);
1461: if (opt) {
1462: TSSetSaveTrajectory(ts);
1463: ts->adjoint_solve = tflg;
1464: }
1465: TSAdjointMonitorSetFromOptions(ts,"-ts_adjoint_monitor","Monitor adjoint timestep size","TSAdjointMonitorDefault",TSAdjointMonitorDefault,NULL);
1466: TSAdjointMonitorSetFromOptions(ts,"-ts_adjoint_monitor_sensi","Monitor sensitivity in the adjoint computation","TSAdjointMonitorSensi",TSAdjointMonitorSensi,NULL);
1467: opt = PETSC_FALSE;
1468: PetscOptionsName("-ts_adjoint_monitor_draw_sensi","Monitor adjoint sensitivities (lambda only) graphically","TSAdjointMonitorDrawSensi",&opt);
1469: if (opt) {
1470: TSMonitorDrawCtx ctx;
1471: PetscInt howoften = 1;
1473: PetscOptionsInt("-ts_adjoint_monitor_draw_sensi","Monitor adjoint sensitivities (lambda only) graphically","TSAdjointMonitorDrawSensi",howoften,&howoften,NULL);
1474: TSMonitorDrawCtxCreate(PetscObjectComm((PetscObject)ts),0,0,PETSC_DECIDE,PETSC_DECIDE,300,300,howoften,&ctx);
1475: TSAdjointMonitorSet(ts,TSAdjointMonitorDrawSensi,ctx,(PetscErrorCode (*)(void**))TSMonitorDrawCtxDestroy);
1476: }
1477: return(0);
1478: }
1480: /*@
1481: TSAdjointStep - Steps one time step backward in the adjoint run
1483: Collective on TS
1485: Input Parameter:
1486: . ts - the TS context obtained from TSCreate()
1488: Level: intermediate
1490: .seealso: TSAdjointSetUp(), TSAdjointSolve()
1491: @*/
1492: PetscErrorCode TSAdjointStep(TS ts)
1493: {
1494: DM dm;
1495: PetscErrorCode ierr;
1499: TSGetDM(ts,&dm);
1500: TSAdjointSetUp(ts);
1502: ts->reason = TS_CONVERGED_ITERATING;
1503: ts->ptime_prev = ts->ptime;
1504: 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);
1505: PetscLogEventBegin(TS_AdjointStep,ts,0,0,0);
1506: (*ts->ops->adjointstep)(ts);
1507: PetscLogEventEnd(TS_AdjointStep,ts,0,0,0);
1508: ts->adjoint_steps++; ts->steps--;
1510: if (ts->reason < 0) {
1511: if (ts->errorifstepfailed) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_NOT_CONVERGED,"TSAdjointStep has failed due to %s",TSConvergedReasons[ts->reason]);
1512: } else if (!ts->reason) {
1513: if (ts->adjoint_steps >= ts->adjoint_max_steps) ts->reason = TS_CONVERGED_ITS;
1514: }
1515: return(0);
1516: }
1518: /*@
1519: TSAdjointSolve - Solves the discrete ajoint problem for an ODE/DAE
1521: Collective on TS
1523: Input Parameter:
1524: . ts - the TS context obtained from TSCreate()
1526: Options Database:
1527: . -ts_adjoint_view_solution <viewerinfo> - views the first gradient with respect to the initial values
1529: Level: intermediate
1531: Notes:
1532: This must be called after a call to TSSolve() that solves the forward problem
1534: By default this will integrate back to the initial time, one can use TSAdjointSetSteps() to step back to a later time
1536: .seealso: TSCreate(), TSSetCostGradients(), TSSetSolution(), TSAdjointStep()
1537: @*/
1538: PetscErrorCode TSAdjointSolve(TS ts)
1539: {
1540: #if defined(TSADJOINT_STAGE)
1541: PetscLogStage adjoint_stage;
1542: #endif
1547: #if defined(TSADJOINT_STAGE)
1548: PetscLogStageRegister("TSAdjoint",&adjoint_stage);
1549: PetscLogStagePush(adjoint_stage);
1550: #endif
1551: TSAdjointSetUp(ts);
1553: /* reset time step and iteration counters */
1554: ts->adjoint_steps = 0;
1555: ts->ksp_its = 0;
1556: ts->snes_its = 0;
1557: ts->num_snes_failures = 0;
1558: ts->reject = 0;
1559: ts->reason = TS_CONVERGED_ITERATING;
1561: if (!ts->adjoint_max_steps) ts->adjoint_max_steps = ts->steps;
1562: if (ts->adjoint_steps >= ts->adjoint_max_steps) ts->reason = TS_CONVERGED_ITS;
1564: while (!ts->reason) {
1565: TSTrajectoryGet(ts->trajectory,ts,ts->steps,&ts->ptime);
1566: TSAdjointMonitor(ts,ts->steps,ts->ptime,ts->vec_sol,ts->numcost,ts->vecs_sensi,ts->vecs_sensip);
1567: TSAdjointEventHandler(ts);
1568: TSAdjointStep(ts);
1569: if ((ts->vec_costintegral || ts->quadraturets) && !ts->costintegralfwd) {
1570: TSAdjointCostIntegral(ts);
1571: }
1572: }
1573: TSTrajectoryGet(ts->trajectory,ts,ts->steps,&ts->ptime);
1574: TSAdjointMonitor(ts,ts->steps,ts->ptime,ts->vec_sol,ts->numcost,ts->vecs_sensi,ts->vecs_sensip);
1575: ts->solvetime = ts->ptime;
1576: TSTrajectoryViewFromOptions(ts->trajectory,NULL,"-ts_trajectory_view");
1577: VecViewFromOptions(ts->vecs_sensi[0],(PetscObject) ts, "-ts_adjoint_view_solution");
1578: ts->adjoint_max_steps = 0;
1579: #if defined(TSADJOINT_STAGE)
1580: PetscLogStagePop();
1581: #endif
1582: return(0);
1583: }
1585: /*@C
1586: TSAdjointMonitor - Runs all user-provided adjoint monitor routines set using TSAdjointMonitorSet()
1588: Collective on TS
1590: Input Parameters:
1591: + ts - time stepping context obtained from TSCreate()
1592: . step - step number that has just completed
1593: . ptime - model time of the state
1594: . u - state at the current model time
1595: . numcost - number of cost functions (dimension of lambda or mu)
1596: . lambda - vectors containing the gradients of the cost functions with respect to the ODE/DAE solution variables
1597: - mu - vectors containing the gradients of the cost functions with respect to the problem parameters
1599: Notes:
1600: TSAdjointMonitor() is typically used automatically within the time stepping implementations.
1601: Users would almost never call this routine directly.
1603: Level: developer
1605: @*/
1606: PetscErrorCode TSAdjointMonitor(TS ts,PetscInt step,PetscReal ptime,Vec u,PetscInt numcost,Vec *lambda, Vec *mu)
1607: {
1609: PetscInt i,n = ts->numberadjointmonitors;
1614: VecLockReadPush(u);
1615: for (i=0; i<n; i++) {
1616: (*ts->adjointmonitor[i])(ts,step,ptime,u,numcost,lambda,mu,ts->adjointmonitorcontext[i]);
1617: }
1618: VecLockReadPop(u);
1619: return(0);
1620: }
1622: /*@
1623: TSAdjointCostIntegral - Evaluate the cost integral in the adjoint run.
1625: Collective on TS
1627: Input Arguments:
1628: . ts - time stepping context
1630: Level: advanced
1632: Notes:
1633: This function cannot be called until TSAdjointStep() has been completed.
1635: .seealso: TSAdjointSolve(), TSAdjointStep
1636: @*/
1637: PetscErrorCode TSAdjointCostIntegral(TS ts)
1638: {
1641: 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);
1642: (*ts->ops->adjointintegral)(ts);
1643: return(0);
1644: }
1646: /* ------------------ Forward (tangent linear) sensitivity ------------------*/
1648: /*@
1649: TSForwardSetUp - Sets up the internal data structures for the later use
1650: of forward sensitivity analysis
1652: Collective on TS
1654: Input Parameter:
1655: . ts - the TS context obtained from TSCreate()
1657: Level: advanced
1659: .seealso: TSCreate(), TSDestroy(), TSSetUp()
1660: @*/
1661: PetscErrorCode TSForwardSetUp(TS ts)
1662: {
1667: if (ts->forwardsetupcalled) return(0);
1668: if (ts->ops->forwardsetup) {
1669: (*ts->ops->forwardsetup)(ts);
1670: }
1671: VecDuplicate(ts->vec_sol,&ts->vec_sensip_col);
1672: ts->forwardsetupcalled = PETSC_TRUE;
1673: return(0);
1674: }
1676: /*@
1677: TSForwardReset - Reset the internal data structures used by forward sensitivity analysis
1679: Collective on TS
1681: Input Parameter:
1682: . ts - the TS context obtained from TSCreate()
1684: Level: advanced
1686: .seealso: TSCreate(), TSDestroy(), TSForwardSetUp()
1687: @*/
1688: PetscErrorCode TSForwardReset(TS ts)
1689: {
1690: TS quadts = ts->quadraturets;
1695: if (ts->ops->forwardreset) {
1696: (*ts->ops->forwardreset)(ts);
1697: }
1698: MatDestroy(&ts->mat_sensip);
1699: if (quadts) {
1700: MatDestroy(&quadts->mat_sensip);
1701: }
1702: VecDestroy(&ts->vec_sensip_col);
1703: ts->forward_solve = PETSC_FALSE;
1704: ts->forwardsetupcalled = PETSC_FALSE;
1705: return(0);
1706: }
1708: /*@
1709: TSForwardSetIntegralGradients - Set the vectors holding forward sensitivities of the integral term.
1711: Input Parameter:
1712: + ts- the TS context obtained from TSCreate()
1713: . numfwdint- number of integrals
1714: - vp = the vectors containing the gradients for each integral w.r.t. parameters
1716: Level: deprecated
1718: .seealso: TSForwardGetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep()
1719: @*/
1720: PetscErrorCode TSForwardSetIntegralGradients(TS ts,PetscInt numfwdint,Vec *vp)
1721: {
1724: 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()");
1725: if (!ts->numcost) ts->numcost = numfwdint;
1727: ts->vecs_integral_sensip = vp;
1728: return(0);
1729: }
1731: /*@
1732: TSForwardGetIntegralGradients - Returns the forward sensitivities ofthe integral term.
1734: Input Parameter:
1735: . ts- the TS context obtained from TSCreate()
1737: Output Parameter:
1738: . vp = the vectors containing the gradients for each integral w.r.t. parameters
1740: Level: deprecated
1742: .seealso: TSForwardSetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep()
1743: @*/
1744: PetscErrorCode TSForwardGetIntegralGradients(TS ts,PetscInt *numfwdint,Vec **vp)
1745: {
1749: if (numfwdint) *numfwdint = ts->numcost;
1750: if (vp) *vp = ts->vecs_integral_sensip;
1751: return(0);
1752: }
1754: /*@
1755: TSForwardStep - Compute the forward sensitivity for one time step.
1757: Collective on TS
1759: Input Arguments:
1760: . ts - time stepping context
1762: Level: advanced
1764: Notes:
1765: This function cannot be called until TSStep() has been completed.
1767: .seealso: TSForwardSetSensitivities(), TSForwardGetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardSetUp()
1768: @*/
1769: PetscErrorCode TSForwardStep(TS ts)
1770: {
1773: if (!ts->ops->forwardstep) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"%s does not provide forward sensitivity analysis",((PetscObject)ts)->type_name);
1774: PetscLogEventBegin(TS_ForwardStep,ts,0,0,0);
1775: (*ts->ops->forwardstep)(ts);
1776: PetscLogEventEnd(TS_ForwardStep,ts,0,0,0);
1777: if (ts->reason < 0 && ts->errorifstepfailed) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_NOT_CONVERGED,"TSFowardStep has failed due to %s",TSConvergedReasons[ts->reason]);
1778: return(0);
1779: }
1781: /*@
1782: TSForwardSetSensitivities - Sets the initial value of the trajectory sensitivities of solution w.r.t. the problem parameters and initial values.
1784: Logically Collective on TS
1786: Input Parameters:
1787: + ts - the TS context obtained from TSCreate()
1788: . nump - number of parameters
1789: - Smat - sensitivities with respect to the parameters, the number of entries in these vectors is the same as the number of parameters
1791: Level: beginner
1793: Notes:
1794: Forward sensitivity is also called 'trajectory sensitivity' in some fields such as power systems.
1795: This function turns on a flag to trigger TSSolve() to compute forward sensitivities automatically.
1796: You must call this function before TSSolve().
1797: The entries in the sensitivity matrix must be correctly initialized with the values S = dy/dp|startingtime.
1799: .seealso: TSForwardGetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep()
1800: @*/
1801: PetscErrorCode TSForwardSetSensitivities(TS ts,PetscInt nump,Mat Smat)
1802: {
1808: ts->forward_solve = PETSC_TRUE;
1809: if (nump == PETSC_DEFAULT) {
1810: MatGetSize(Smat,NULL,&ts->num_parameters);
1811: } else ts->num_parameters = nump;
1812: PetscObjectReference((PetscObject)Smat);
1813: MatDestroy(&ts->mat_sensip);
1814: ts->mat_sensip = Smat;
1815: return(0);
1816: }
1818: /*@
1819: TSForwardGetSensitivities - Returns the trajectory sensitivities
1821: Not Collective, but Vec returned is parallel if TS is parallel
1823: Output Parameter:
1824: + ts - the TS context obtained from TSCreate()
1825: . nump - number of parameters
1826: - Smat - sensitivities with respect to the parameters, the number of entries in these vectors is the same as the number of parameters
1828: Level: intermediate
1830: .seealso: TSForwardSetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep()
1831: @*/
1832: PetscErrorCode TSForwardGetSensitivities(TS ts,PetscInt *nump,Mat *Smat)
1833: {
1836: if (nump) *nump = ts->num_parameters;
1837: if (Smat) *Smat = ts->mat_sensip;
1838: return(0);
1839: }
1841: /*@
1842: TSForwardCostIntegral - Evaluate the cost integral in the forward run.
1844: Collective on TS
1846: Input Arguments:
1847: . ts - time stepping context
1849: Level: advanced
1851: Notes:
1852: This function cannot be called until TSStep() has been completed.
1854: .seealso: TSSolve(), TSAdjointCostIntegral()
1855: @*/
1856: PetscErrorCode TSForwardCostIntegral(TS ts)
1857: {
1860: 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);
1861: (*ts->ops->forwardintegral)(ts);
1862: return(0);
1863: }
1865: /*@
1866: TSForwardSetInitialSensitivities - Set initial values for tangent linear sensitivities
1868: Collective on TS
1870: Input Parameter
1871: + ts - the TS context obtained from TSCreate()
1872: - didp - parametric sensitivities of the initial condition
1874: Level: intermediate
1876: 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.
1878: .seealso: TSForwardSetSensitivities()
1879: @*/
1880: PetscErrorCode TSForwardSetInitialSensitivities(TS ts,Mat didp)
1881: {
1886: if (!ts->mat_sensip) {
1887: TSForwardSetSensitivities(ts,PETSC_DEFAULT,didp);
1888: }
1889: return(0);
1890: }
1892: /*@
1893: TSForwardGetStages - Get the number of stages and the tangent linear sensitivities at the intermediate stages
1895: Input Parameter:
1896: . ts - the TS context obtained from TSCreate()
1898: Output Parameters:
1899: + ns - number of stages
1900: - S - tangent linear sensitivities at the intermediate stages
1902: Level: advanced
1904: @*/
1905: PetscErrorCode TSForwardGetStages(TS ts,PetscInt *ns,Mat **S)
1906: {
1912: if (!ts->ops->getstages) *S=NULL;
1913: else {
1914: (*ts->ops->forwardgetstages)(ts,ns,S);
1915: }
1916: return(0);
1917: }
1919: /*@
1920: TSCreateQuadratureTS - Create a sub-TS that evaluates integrals over time
1922: Input Parameter:
1923: + ts - the TS context obtained from TSCreate()
1924: - fwd - flag indicating whether to evaluate cost integral in the forward run or the adjoint run
1926: Output Parameters:
1927: . quadts - the child TS context
1929: Level: intermediate
1931: .seealso: TSGetQuadratureTS()
1932: @*/
1933: PetscErrorCode TSCreateQuadratureTS(TS ts,PetscBool fwd,TS *quadts)
1934: {
1935: char prefix[128];
1941: TSDestroy(&ts->quadraturets);
1942: TSCreate(PetscObjectComm((PetscObject)ts),&ts->quadraturets);
1943: PetscObjectIncrementTabLevel((PetscObject)ts->quadraturets,(PetscObject)ts,1);
1944: PetscLogObjectParent((PetscObject)ts,(PetscObject)ts->quadraturets);
1945: PetscSNPrintf(prefix,sizeof(prefix),"%squad_",((PetscObject)ts)->prefix ? ((PetscObject)ts)->prefix : "");
1946: TSSetOptionsPrefix(ts->quadraturets,prefix);
1947: *quadts = ts->quadraturets;
1949: if (ts->numcost) {
1950: VecCreateSeq(PETSC_COMM_SELF,ts->numcost,&(*quadts)->vec_sol);
1951: } else {
1952: VecCreateSeq(PETSC_COMM_SELF,1,&(*quadts)->vec_sol);
1953: }
1954: ts->costintegralfwd = fwd;
1955: return(0);
1956: }
1958: /*@
1959: TSGetQuadratureTS - Return the sub-TS that evaluates integrals over time
1961: Input Parameter:
1962: . ts - the TS context obtained from TSCreate()
1964: Output Parameters:
1965: + fwd - flag indicating whether to evaluate cost integral in the forward run or the adjoint run
1966: - quadts - the child TS context
1968: Level: intermediate
1970: .seealso: TSCreateQuadratureTS()
1971: @*/
1972: PetscErrorCode TSGetQuadratureTS(TS ts,PetscBool *fwd,TS *quadts)
1973: {
1976: if (fwd) *fwd = ts->costintegralfwd;
1977: if (quadts) *quadts = ts->quadraturets;
1978: return(0);
1979: }
1981: /*@
1982: TSComputeSNESJacobian - Compute the SNESJacobian
1984: Input Parameters:
1985: + ts - the TS context obtained from TSCreate()
1986: - x - state vector
1988: Output Parameters:
1989: + J - Jacobian matrix
1990: - Jpre - preconditioning matrix for J (may be same as J)
1992: Level: developer
1994: Notes:
1995: Using SNES to compute the Jacobian enables finite differencing when TS Jacobian is not available.
1996: @*/
1997: PetscErrorCode TSComputeSNESJacobian(TS ts,Vec x,Mat J,Mat Jpre)
1998: {
1999: SNES snes = ts->snes;
2000: PetscErrorCode (*jac)(SNES,Vec,Mat,Mat,void*) = NULL;
2004: /*
2005: Unlike implicit methods, explicit methods do not have SNESMatFDColoring in the snes object
2006: because SNESSolve() has not been called yet; so querying SNESMatFDColoring does not work for
2007: explicit methods. Instead, we check the Jacobian compute function directly to determin if FD
2008: coloring is used.
2009: */
2010: SNESGetJacobian(snes,NULL,NULL,&jac,NULL);
2011: if (jac == SNESComputeJacobianDefaultColor) {
2012: Vec f;
2013: SNESSetSolution(snes,x);
2014: SNESGetFunction(snes,&f,NULL,NULL);
2015: /* Force MatFDColoringApply to evaluate the SNES residual function for the base vector */
2016: SNESComputeFunction(snes,x,f);
2017: }
2018: SNESComputeJacobian(snes,x,J,Jpre);
2019: return(0);
2020: }