Actual source code: tsimpl.h
1: #ifndef __TSIMPL_H
4: #include <petscts.h>
5: #include <petsc/private/petscimpl.h>
7: /*
8: Timesteping context.
9: General DAE: F(t,U,U_t) = 0, required Jacobian is G'(U) where G(U) = F(t,U,U0+a*U)
10: General ODE: U_t = F(t,U) <-- the right-hand-side function
11: Linear ODE: U_t = A(t) U <-- the right-hand-side matrix
12: Linear (no time) ODE: U_t = A U <-- the right-hand-side matrix
13: */
15: /*
16: Maximum number of monitors you can run with a single TS
17: */
18: #define MAXTSMONITORS 10
20: PETSC_EXTERN PetscBool TSRegisterAllCalled;
21: PETSC_EXTERN PetscErrorCode TSRegisterAll(void);
22: PETSC_EXTERN PetscErrorCode TSAdaptRegisterAll(void);
24: PETSC_EXTERN PetscErrorCode TSRKRegisterAll(void);
25: PETSC_EXTERN PetscErrorCode TSMPRKRegisterAll(void);
26: PETSC_EXTERN PetscErrorCode TSARKIMEXRegisterAll(void);
27: PETSC_EXTERN PetscErrorCode TSRosWRegisterAll(void);
28: PETSC_EXTERN PetscErrorCode TSGLLERegisterAll(void);
29: PETSC_EXTERN PetscErrorCode TSGLLEAdaptRegisterAll(void);
30: PETSC_EXTERN PetscErrorCode TSIRKRegisterAll(void);
32: typedef struct _TSOps *TSOps;
34: struct _TSOps {
35: PetscErrorCode (*snesfunction)(SNES, Vec, Vec, TS);
36: PetscErrorCode (*snesjacobian)(SNES, Vec, Mat, Mat, TS);
37: PetscErrorCode (*setup)(TS);
38: PetscErrorCode (*step)(TS);
39: PetscErrorCode (*solve)(TS);
40: PetscErrorCode (*interpolate)(TS, PetscReal, Vec);
41: PetscErrorCode (*evaluatewlte)(TS, NormType, PetscInt *, PetscReal *);
42: PetscErrorCode (*evaluatestep)(TS, PetscInt, Vec, PetscBool *);
43: PetscErrorCode (*setfromoptions)(TS, PetscOptionItems *);
44: PetscErrorCode (*destroy)(TS);
45: PetscErrorCode (*view)(TS, PetscViewer);
46: PetscErrorCode (*reset)(TS);
47: PetscErrorCode (*linearstability)(TS, PetscReal, PetscReal, PetscReal *, PetscReal *);
48: PetscErrorCode (*load)(TS, PetscViewer);
49: PetscErrorCode (*rollback)(TS);
50: PetscErrorCode (*getstages)(TS, PetscInt *, Vec *[]);
51: PetscErrorCode (*adjointstep)(TS);
52: PetscErrorCode (*adjointsetup)(TS);
53: PetscErrorCode (*adjointreset)(TS);
54: PetscErrorCode (*adjointintegral)(TS);
55: PetscErrorCode (*forwardsetup)(TS);
56: PetscErrorCode (*forwardreset)(TS);
57: PetscErrorCode (*forwardstep)(TS);
58: PetscErrorCode (*forwardintegral)(TS);
59: PetscErrorCode (*forwardgetstages)(TS, PetscInt *, Mat *[]);
60: PetscErrorCode (*getsolutioncomponents)(TS, PetscInt *, Vec *);
61: PetscErrorCode (*getauxsolution)(TS, Vec *);
62: PetscErrorCode (*gettimeerror)(TS, PetscInt, Vec *);
63: PetscErrorCode (*settimeerror)(TS, Vec);
64: PetscErrorCode (*startingmethod)(TS);
65: PetscErrorCode (*initcondition)(TS, Vec);
66: PetscErrorCode (*exacterror)(TS, Vec, Vec);
67: };
69: /*
70: TSEvent - Abstract object to handle event monitoring
71: */
72: typedef struct _n_TSEvent *TSEvent;
74: typedef struct _TSTrajectoryOps *TSTrajectoryOps;
76: struct _TSTrajectoryOps {
77: PetscErrorCode (*view)(TSTrajectory, PetscViewer);
78: PetscErrorCode (*reset)(TSTrajectory);
79: PetscErrorCode (*destroy)(TSTrajectory);
80: PetscErrorCode (*set)(TSTrajectory, TS, PetscInt, PetscReal, Vec);
81: PetscErrorCode (*get)(TSTrajectory, TS, PetscInt, PetscReal *);
82: PetscErrorCode (*setfromoptions)(TSTrajectory, PetscOptionItems *);
83: PetscErrorCode (*setup)(TSTrajectory, TS);
84: };
86: /* TSHistory is an helper object that allows inquiring
87: the TSTrajectory by time and not by the step number only */
88: typedef struct _n_TSHistory *TSHistory;
90: struct _p_TSTrajectory {
91: PETSCHEADER(struct _TSTrajectoryOps);
92: TSHistory tsh; /* associates times to unique step ids */
93: /* stores necessary data to reconstruct states and derivatives via Lagrangian interpolation */
94: struct {
95: PetscInt order; /* interpolation order */
96: Vec *W; /* work vectors */
97: PetscScalar *L; /* workspace for Lagrange basis */
98: PetscReal *T; /* Lagrange times (stored) */
99: Vec *WW; /* just an array of pointers */
100: PetscBool *TT; /* workspace for Lagrange */
101: PetscReal *TW; /* Lagrange times (workspace) */
103: /* caching */
104: PetscBool caching;
105: struct {
106: PetscObjectId id;
107: PetscObjectState state;
108: PetscReal time;
109: PetscInt step;
110: } Ucached;
111: struct {
112: PetscObjectId id;
113: PetscObjectState state;
114: PetscReal time;
115: PetscInt step;
116: } Udotcached;
117: } lag;
118: Vec U, Udot; /* used by TSTrajectory{Get|Restore}UpdatedHistoryVecs */
119: PetscBool usehistory; /* whether to use TSHistory */
120: PetscBool solution_only; /* whether we dump just the solution or also the stages */
121: PetscBool adjoint_solve_mode; /* whether we will use the Trajectory inside a TSAdjointSolve() or not */
122: PetscViewer monitor;
123: PetscInt setupcalled; /* true if setup has been called */
124: PetscInt recomps; /* counter for recomputations in the adjoint run */
125: PetscInt diskreads, diskwrites; /* counters for disk checkpoint reads and writes */
126: char **names; /* the name of each variable; each process has only the local names */
127: PetscBool keepfiles; /* keep the files generated during the run after the run is complete */
128: char *dirname, *filetemplate; /* directory name and file name template for disk checkpoints */
129: char *dirfiletemplate; /* complete directory and file name template for disk checkpoints */
130: PetscErrorCode (*transform)(void *, Vec, Vec *);
131: PetscErrorCode (*transformdestroy)(void *);
132: void *transformctx;
133: void *data;
134: };
136: typedef struct _TS_RHSSplitLink *TS_RHSSplitLink;
137: struct _TS_RHSSplitLink {
138: TS ts;
139: char *splitname;
140: IS is;
141: TS_RHSSplitLink next;
142: PetscLogEvent event;
143: };
145: typedef struct _TS_TimeSpan *TSTimeSpan;
146: struct _TS_TimeSpan {
147: PetscInt num_span_times; /* number of time points */
148: PetscReal *span_times; /* array of the time span */
149: PetscReal reltol; /* relative tolerance for span point detection */
150: PetscReal abstol; /* absolute tolerance for span point detection */
151: PetscInt spanctr; /* counter of the time points that have been reached */
152: Vec *vecs_sol; /* array of the solutions at the specified time points */
153: };
155: struct _p_TS {
156: PETSCHEADER(struct _TSOps);
157: TSProblemType problem_type;
158: TSEquationType equation_type;
160: DM dm;
161: Vec vec_sol; /* solution vector in first and second order equations */
162: Vec vec_dot; /* time derivative vector in second order equations */
163: TSAdapt adapt;
164: TSAdaptType default_adapt_type;
165: TSEvent event;
167: /* ---------------- User (or PETSc) Provided stuff ---------------------*/
168: PetscErrorCode (*monitor[MAXTSMONITORS])(TS, PetscInt, PetscReal, Vec, void *);
169: PetscErrorCode (*monitordestroy[MAXTSMONITORS])(void **);
170: void *monitorcontext[MAXTSMONITORS];
171: PetscInt numbermonitors;
172: PetscErrorCode (*adjointmonitor[MAXTSMONITORS])(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *, void *);
173: PetscErrorCode (*adjointmonitordestroy[MAXTSMONITORS])(void **);
174: void *adjointmonitorcontext[MAXTSMONITORS];
175: PetscInt numberadjointmonitors;
176: PetscInt monitorFrequency; /* Number of timesteps between monitor output */
178: PetscErrorCode (*prestep)(TS);
179: PetscErrorCode (*prestage)(TS, PetscReal);
180: PetscErrorCode (*poststage)(TS, PetscReal, PetscInt, Vec *);
181: PetscErrorCode (*postevaluate)(TS);
182: PetscErrorCode (*poststep)(TS);
183: PetscErrorCode (*functiondomainerror)(TS, PetscReal, Vec, PetscBool *);
185: /* ---------------------- Sensitivity Analysis support -----------------*/
186: TSTrajectory trajectory; /* All solutions are kept here for the entire time integration process */
187: Vec *vecs_sensi; /* one vector for each cost function */
188: Vec *vecs_sensip;
189: PetscInt numcost; /* number of cost functions */
190: Vec vec_costintegral;
191: PetscInt adjointsetupcalled;
192: PetscInt adjoint_steps;
193: PetscInt adjoint_max_steps;
194: PetscBool adjoint_solve; /* immediately call TSAdjointSolve() after TSSolve() is complete */
195: PetscBool costintegralfwd; /* cost integral is evaluated in the forward run if true */
196: Vec vec_costintegrand; /* workspace for Adjoint computations */
197: Mat Jacp, Jacprhs;
198: void *ijacobianpctx, *rhsjacobianpctx;
199: void *costintegrandctx;
200: Vec *vecs_drdu;
201: Vec *vecs_drdp;
202: Vec vec_drdu_col, vec_drdp_col;
204: /* first-order adjoint */
205: PetscErrorCode (*rhsjacobianp)(TS, PetscReal, Vec, Mat, void *);
206: PetscErrorCode (*ijacobianp)(TS, PetscReal, Vec, Vec, PetscReal, Mat, void *);
207: PetscErrorCode (*costintegrand)(TS, PetscReal, Vec, Vec, void *);
208: PetscErrorCode (*drdufunction)(TS, PetscReal, Vec, Vec *, void *);
209: PetscErrorCode (*drdpfunction)(TS, PetscReal, Vec, Vec *, void *);
211: /* second-order adjoint */
212: Vec *vecs_sensi2;
213: Vec *vecs_sensi2p;
214: Vec vec_dir; /* directional vector for optimization */
215: Vec *vecs_fuu, *vecs_guu;
216: Vec *vecs_fup, *vecs_gup;
217: Vec *vecs_fpu, *vecs_gpu;
218: Vec *vecs_fpp, *vecs_gpp;
219: void *ihessianproductctx, *rhshessianproductctx;
220: PetscErrorCode (*ihessianproduct_fuu)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *);
221: PetscErrorCode (*ihessianproduct_fup)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *);
222: PetscErrorCode (*ihessianproduct_fpu)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *);
223: PetscErrorCode (*ihessianproduct_fpp)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *);
224: PetscErrorCode (*rhshessianproduct_guu)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *);
225: PetscErrorCode (*rhshessianproduct_gup)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *);
226: PetscErrorCode (*rhshessianproduct_gpu)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *);
227: PetscErrorCode (*rhshessianproduct_gpp)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *);
229: /* specific to forward sensitivity analysis */
230: Mat mat_sensip; /* matrix storing forward sensitivities */
231: Vec vec_sensip_col; /* space for a column of the sensip matrix */
232: Vec *vecs_integral_sensip; /* one vector for each integral */
233: PetscInt num_parameters;
234: PetscInt num_initialvalues;
235: void *vecsrhsjacobianpctx;
236: PetscInt forwardsetupcalled;
237: PetscBool forward_solve;
238: PetscErrorCode (*vecsrhsjacobianp)(TS, PetscReal, Vec, Vec *, void *);
240: /* ---------------------- IMEX support ---------------------------------*/
241: /* These extra slots are only used when the user provides both Implicit and RHS */
242: Mat Arhs; /* Right hand side matrix */
243: Mat Brhs; /* Right hand side preconditioning matrix */
244: Vec Frhs; /* Right hand side function value */
246: /* This is a general caching scheme to avoid recomputing the Jacobian at a place that has been previously been evaluated.
247: * The present use case is that TSComputeRHSFunctionLinear() evaluates the Jacobian once and we don't want it to be immeditely re-evaluated.
248: */
249: struct {
250: PetscReal time; /* The time at which the matrices were last evaluated */
251: PetscObjectId Xid; /* Unique ID of solution vector at which the Jacobian was last evaluated */
252: PetscObjectState Xstate; /* State of the solution vector */
253: MatStructure mstructure; /* The structure returned */
254: /* Flag to unshift Jacobian before calling the IJacobian or RHSJacobian functions. This is useful
255: * if the user would like to reuse (part of) the Jacobian from the last evaluation. */
256: PetscBool reuse;
257: PetscReal scale, shift;
258: } rhsjacobian;
260: struct {
261: PetscReal shift; /* The derivative of the lhs wrt to Xdot */
262: } ijacobian;
264: MatStructure axpy_pattern; /* information about the nonzero pattern of the RHS Jacobian in reference to the implicit Jacobian */
265: /* --------------------Nonlinear Iteration------------------------------*/
266: SNES snes;
267: PetscBool usessnes; /* Flag set by each TSType to indicate if the type actually uses a SNES;
268: this works around the design flaw that a SNES is ALWAYS created with TS even when it is not needed.*/
269: PetscInt ksp_its; /* total number of linear solver iterations */
270: PetscInt snes_its; /* total number of nonlinear solver iterations */
271: PetscInt num_snes_failures;
272: PetscInt max_snes_failures;
274: /* --- Logging --- */
275: PetscInt ifuncs, rhsfuncs, ijacs, rhsjacs;
277: /* --- Data that is unique to each particular solver --- */
278: PetscInt setupcalled; /* true if setup has been called */
279: void *data; /* implementationspecific data */
280: void *user; /* user context */
282: /* ------------------ Parameters -------------------------------------- */
283: PetscInt max_steps; /* max number of steps */
284: PetscReal max_time; /* max time allowed */
286: /* --------------------------------------------------------------------- */
288: PetscBool steprollback; /* flag to indicate that the step was rolled back */
289: PetscBool steprestart; /* flag to indicate that the timestepper has to discard any history and restart */
290: PetscInt steps; /* steps taken so far in all successive calls to TSSolve() */
291: PetscReal ptime; /* time at the start of the current step (stage time is internal if it exists) */
292: PetscReal time_step; /* current time increment */
293: PetscReal ptime_prev; /* time at the start of the previous step */
294: PetscReal ptime_prev_rollback; /* time at the start of the 2nd previous step to recover from rollback */
295: PetscReal solvetime; /* time at the conclusion of TSSolve() */
296: PetscBool stifflyaccurate; /* flag to indicate that the method is stiffly accurate */
298: TSConvergedReason reason;
299: PetscBool errorifstepfailed;
300: PetscInt reject, max_reject;
301: TSExactFinalTimeOption exact_final_time;
303: PetscReal atol, rtol; /* Relative and absolute tolerance for local truncation error */
304: Vec vatol, vrtol; /* Relative and absolute tolerance in vector form */
305: PetscReal cfltime, cfltime_local;
307: PetscBool testjacobian;
308: PetscBool testjacobiantranspose;
309: /* ------------------- Default work-area management ------------------ */
310: PetscInt nwork;
311: Vec *work;
313: /* ---------------------- RHS splitting support ---------------------------------*/
314: PetscInt num_rhs_splits;
315: TS_RHSSplitLink tsrhssplit;
316: PetscBool use_splitrhsfunction;
318: /* ---------------------- Quadrature integration support ---------------------------------*/
319: TS quadraturets;
321: /* ---------------------- Time span support ---------------------------------*/
322: TSTimeSpan tspan;
323: };
325: struct _TSAdaptOps {
326: PetscErrorCode (*choose)(TSAdapt, TS, PetscReal, PetscInt *, PetscReal *, PetscBool *, PetscReal *, PetscReal *, PetscReal *);
327: PetscErrorCode (*destroy)(TSAdapt);
328: PetscErrorCode (*reset)(TSAdapt);
329: PetscErrorCode (*view)(TSAdapt, PetscViewer);
330: PetscErrorCode (*setfromoptions)(TSAdapt, PetscOptionItems *);
331: PetscErrorCode (*load)(TSAdapt, PetscViewer);
332: };
334: struct _p_TSAdapt {
335: PETSCHEADER(struct _TSAdaptOps);
336: void *data;
337: PetscErrorCode (*checkstage)(TSAdapt, TS, PetscReal, Vec, PetscBool *);
338: struct {
339: PetscInt n; /* number of candidate schemes, including the one currently in use */
340: PetscBool inuse_set; /* the current scheme has been set */
341: const char *name[16]; /* name of the scheme */
342: PetscInt order[16]; /* classical order of each scheme */
343: PetscInt stageorder[16]; /* stage order of each scheme */
344: PetscReal ccfl[16]; /* stability limit relative to explicit Euler */
345: PetscReal cost[16]; /* relative measure of the amount of work required for each scheme */
346: } candidates;
347: PetscBool always_accept;
348: PetscReal safety; /* safety factor relative to target error/stability goal */
349: PetscReal reject_safety; /* extra safety factor if the last step was rejected */
350: PetscReal clip[2]; /* admissible time step decrease/increase factors */
351: PetscReal dt_min, dt_max; /* admissible minimum and maximum time step */
352: PetscReal ignore_max; /* minimum value of the solution to be considered by the adaptor */
353: PetscBool glee_use_local; /* GLEE adaptor uses global or local error */
354: PetscReal scale_solve_failed; /* scale step by this factor if solver (linear or nonlinear) fails. */
355: PetscReal matchstepfac[2]; /* factors to control the behaviour of matchstep */
356: NormType wnormtype;
357: PetscViewer monitor;
358: PetscInt timestepjustdecreased_delay; /* number of timesteps after a decrease in the timestep before the timestep can be increased */
359: PetscInt timestepjustdecreased;
360: PetscReal dt_span_cached; /* time step before hitting a TS span time point */
361: };
363: typedef struct _p_DMTS *DMTS;
364: typedef struct _DMTSOps *DMTSOps;
365: struct _DMTSOps {
366: TSRHSFunction rhsfunction;
367: TSRHSJacobian rhsjacobian;
369: TSIFunction ifunction;
370: PetscErrorCode (*ifunctionview)(void *, PetscViewer);
371: PetscErrorCode (*ifunctionload)(void **, PetscViewer);
373: TSIJacobian ijacobian;
374: PetscErrorCode (*ijacobianview)(void *, PetscViewer);
375: PetscErrorCode (*ijacobianload)(void **, PetscViewer);
377: TSI2Function i2function;
378: TSI2Jacobian i2jacobian;
380: TSTransientVariable transientvar;
382: TSSolutionFunction solution;
383: TSForcingFunction forcing;
385: PetscErrorCode (*destroy)(DMTS);
386: PetscErrorCode (*duplicate)(DMTS, DMTS);
387: };
389: struct _p_DMTS {
390: PETSCHEADER(struct _DMTSOps);
391: PetscContainer rhsfunctionctxcontainer;
392: PetscContainer rhsjacobianctxcontainer;
394: PetscContainer ifunctionctxcontainer;
395: PetscContainer ijacobianctxcontainer;
397: PetscContainer i2functionctxcontainer;
398: PetscContainer i2jacobianctxcontainer;
400: void *transientvarctx;
402: void *solutionctx;
403: void *forcingctx;
405: void *data;
407: /* This is NOT reference counted. The DM on which this context was first created is cached here to implement one-way
408: * copy-on-write. When DMGetDMTSWrite() sees a request using a different DM, it makes a copy. Thus, if a user
409: * only interacts directly with one level, e.g., using TSSetIFunction(), then coarse levels of a multilevel item
410: * integrator are built, then the user changes the routine with another call to TSSetIFunction(), it automatically
411: * propagates to all the levels. If instead, they get out a specific level and set the function on that level,
412: * subsequent changes to the original level will no longer propagate to that level.
413: */
414: DM originaldm;
415: };
417: PETSC_EXTERN PetscErrorCode DMTSUnsetRHSFunctionContext_Internal(DM);
418: PETSC_EXTERN PetscErrorCode DMTSUnsetRHSJacobianContext_Internal(DM);
419: PETSC_EXTERN PetscErrorCode DMTSUnsetIFunctionContext_Internal(DM);
420: PETSC_EXTERN PetscErrorCode DMTSUnsetIJacobianContext_Internal(DM);
421: PETSC_EXTERN PetscErrorCode DMTSUnsetI2FunctionContext_Internal(DM);
422: PETSC_EXTERN PetscErrorCode DMTSUnsetI2JacobianContext_Internal(DM);
424: PETSC_EXTERN PetscErrorCode DMGetDMTS(DM, DMTS *);
425: PETSC_EXTERN PetscErrorCode DMGetDMTSWrite(DM, DMTS *);
426: PETSC_EXTERN PetscErrorCode DMCopyDMTS(DM, DM);
427: PETSC_EXTERN PetscErrorCode DMTSView(DMTS, PetscViewer);
428: PETSC_EXTERN PetscErrorCode DMTSLoad(DMTS, PetscViewer);
429: PETSC_EXTERN PetscErrorCode DMTSCopy(DMTS, DMTS);
431: typedef enum {
432: TSEVENT_NONE,
433: TSEVENT_LOCATED_INTERVAL,
434: TSEVENT_PROCESSING,
435: TSEVENT_ZERO,
436: TSEVENT_RESET_NEXTSTEP
437: } TSEventStatus;
439: struct _n_TSEvent {
440: PetscScalar *fvalue; /* value of event function at the end of the step*/
441: PetscScalar *fvalue_prev; /* value of event function at start of the step (left end-point of event interval) */
442: PetscReal ptime_prev; /* time at step start (left end-point of event interval) */
443: PetscReal ptime_end; /* end time of step (when an event interval is detected, ptime_end is fixed to the time at step end during event processing) */
444: PetscReal ptime_right; /* time on the right end-point of the event interval */
445: PetscScalar *fvalue_right; /* value of event function at the right end-point of the event interval */
446: PetscInt *side; /* Used for detecting repetition of end-point, -1 => left, +1 => right */
447: PetscReal timestep_prev; /* previous time step */
448: PetscReal timestep_posteventinterval; /* time step immediately after the event interval */
449: PetscReal timestep_postevent; /* time step immediately after the event */
450: PetscReal timestep_min; /* Minimum time step */
451: PetscBool *zerocrossing; /* Flag to signal zero crossing detection */
452: PetscErrorCode (*eventhandler)(TS, PetscReal, Vec, PetscScalar *, void *); /* User event handler function */
453: PetscErrorCode (*postevent)(TS, PetscInt, PetscInt[], PetscReal, Vec, PetscBool, void *); /* User post event function */
454: void *ctx; /* User context for event handler and post even functions */
455: PetscInt *direction; /* Zero crossing direction: 1 -> Going positive, -1 -> Going negative, 0 -> Any */
456: PetscBool *terminate; /* 1 -> Terminate time stepping, 0 -> continue */
457: PetscInt nevents; /* Number of events to handle */
458: PetscInt nevents_zero; /* Number of event zero detected */
459: PetscInt *events_zero; /* List of events that have reached zero */
460: PetscReal *vtol; /* Vector tolerances for event zero check */
461: TSEventStatus status; /* Event status */
462: PetscInt iterctr; /* Iteration counter */
463: PetscViewer monitor;
464: /* Struct to record the events */
465: struct {
466: PetscInt ctr; /* recorder counter */
467: PetscReal *time; /* Event times */
468: PetscInt *stepnum; /* Step numbers */
469: PetscInt *nevents; /* Number of events occurring at the event times */
470: PetscInt **eventidx; /* Local indices of the events in the event list */
471: } recorder;
472: PetscInt recsize; /* Size of recorder stack */
473: PetscInt refct; /* reference count */
474: };
476: PETSC_EXTERN PetscErrorCode TSEventInitialize(TSEvent, TS, PetscReal, Vec);
477: PETSC_EXTERN PetscErrorCode TSEventDestroy(TSEvent *);
478: PETSC_EXTERN PetscErrorCode TSEventHandler(TS);
479: PETSC_EXTERN PetscErrorCode TSAdjointEventHandler(TS);
481: PETSC_EXTERN PetscLogEvent TS_AdjointStep;
482: PETSC_EXTERN PetscLogEvent TS_Step;
483: PETSC_EXTERN PetscLogEvent TS_PseudoComputeTimeStep;
484: PETSC_EXTERN PetscLogEvent TS_FunctionEval;
485: PETSC_EXTERN PetscLogEvent TS_JacobianEval;
486: PETSC_EXTERN PetscLogEvent TS_ForwardStep;
488: typedef enum {
489: TS_STEP_INCOMPLETE, /* vec_sol, ptime, etc point to beginning of step */
490: TS_STEP_PENDING, /* vec_sol advanced, but step has not been accepted yet */
491: TS_STEP_COMPLETE /* step accepted and ptime, steps, etc have been advanced */
492: } TSStepStatus;
494: struct _n_TSMonitorLGCtx {
495: PetscDrawLG lg;
496: PetscBool semilogy;
497: PetscInt howoften; /* when > 0 uses step % howoften, when negative only final solution plotted */
498: PetscInt ksp_its, snes_its;
499: char **names;
500: char **displaynames;
501: PetscInt ndisplayvariables;
502: PetscInt *displayvariables;
503: PetscReal *displayvalues;
504: PetscErrorCode (*transform)(void *, Vec, Vec *);
505: PetscErrorCode (*transformdestroy)(void *);
506: void *transformctx;
507: };
509: struct _n_TSMonitorSPCtx {
510: PetscDrawSP sp;
511: PetscInt howoften; /* when > 0 uses step % howoften, when negative only final solution plotted */
512: PetscInt retain; /* Retain n points plotted to show trajectories, or -1 for all points */
513: PetscBool phase; /* Plot in phase space rather than coordinate space */
514: PetscInt ksp_its, snes_its;
515: };
517: struct _n_TSMonitorEnvelopeCtx {
518: Vec max, min;
519: };
521: /*
522: Checks if the user provide a TSSetIFunction() but an explicit method is called; generate an error in that case
523: */
524: static inline PetscErrorCode TSCheckImplicitTerm(TS ts)
525: {
526: TSIFunction ifunction;
527: DM dm;
529: TSGetDM(ts, &dm);
530: DMTSGetIFunction(dm, &ifunction, NULL);
532: return 0;
533: }
535: PETSC_EXTERN PetscErrorCode TSGetRHSMats_Private(TS, Mat *, Mat *);
536: /* this is declared here as TSHistory is not public */
537: PETSC_EXTERN PetscErrorCode TSAdaptHistorySetTSHistory(TSAdapt, TSHistory, PetscBool);
539: PETSC_INTERN PetscErrorCode TSTrajectoryReconstruct_Private(TSTrajectory, TS, PetscReal, Vec, Vec);
540: PETSC_INTERN PetscErrorCode TSTrajectorySetUp_Basic(TSTrajectory, TS);
542: PETSC_EXTERN PetscLogEvent TSTrajectory_Set;
543: PETSC_EXTERN PetscLogEvent TSTrajectory_Get;
544: PETSC_EXTERN PetscLogEvent TSTrajectory_GetVecs;
545: PETSC_EXTERN PetscLogEvent TSTrajectory_SetUp;
546: PETSC_EXTERN PetscLogEvent TSTrajectory_DiskWrite;
547: PETSC_EXTERN PetscLogEvent TSTrajectory_DiskRead;
549: struct _n_TSMonitorDrawCtx {
550: PetscViewer viewer;
551: Vec initialsolution;
552: PetscBool showinitial;
553: PetscInt howoften; /* when > 0 uses step % howoften, when negative only final solution plotted */
554: PetscBool showtimestepandtime;
555: };
556: #endif