Actual source code: tsimpl.h
1: #pragma once
3: #include <petscts.h>
4: #include <petsc/private/petscimpl.h>
6: /* SUBMANSEC = TS */
8: /*
9: Timesteping context.
10: General DAE: F(t,U,U_t) = 0, required Jacobian is G'(U) where G(U) = F(t,U,U0+a*U)
11: General ODE: U_t = F(t,U) <-- the right-hand-side function
12: Linear ODE: U_t = A(t) U <-- the right-hand-side matrix
13: Linear (no time) ODE: U_t = A U <-- the right-hand-side matrix
14: */
16: /*
17: Maximum number of monitors you can run with a single TS
18: */
19: #define MAXTSMONITORS 10
21: PETSC_EXTERN PetscBool TSRegisterAllCalled;
22: PETSC_EXTERN PetscErrorCode TSRegisterAll(void);
23: PETSC_EXTERN PetscErrorCode TSAdaptRegisterAll(void);
25: PETSC_EXTERN PetscErrorCode TSRKRegisterAll(void);
26: PETSC_EXTERN PetscErrorCode TSMPRKRegisterAll(void);
27: PETSC_EXTERN PetscErrorCode TSARKIMEXRegisterAll(void);
28: PETSC_EXTERN PetscErrorCode TSRosWRegisterAll(void);
29: PETSC_EXTERN PetscErrorCode TSGLLERegisterAll(void);
30: PETSC_EXTERN PetscErrorCode TSGLLEAdaptRegisterAll(void);
31: PETSC_EXTERN PetscErrorCode TSIRKRegisterAll(void);
33: typedef struct _TSOps *TSOps;
35: struct _TSOps {
36: PetscErrorCode (*snesfunction)(SNES, Vec, Vec, TS);
37: PetscErrorCode (*snesjacobian)(SNES, Vec, Mat, Mat, TS);
38: PetscErrorCode (*setup)(TS);
39: PetscErrorCode (*step)(TS);
40: PetscErrorCode (*solve)(TS);
41: PetscErrorCode (*interpolate)(TS, PetscReal, Vec);
42: PetscErrorCode (*evaluatewlte)(TS, NormType, PetscInt *, PetscReal *);
43: PetscErrorCode (*evaluatestep)(TS, PetscInt, Vec, PetscBool *);
44: PetscErrorCode (*setfromoptions)(TS, PetscOptionItems *);
45: PetscErrorCode (*destroy)(TS);
46: PetscErrorCode (*view)(TS, PetscViewer);
47: PetscErrorCode (*reset)(TS);
48: PetscErrorCode (*linearstability)(TS, PetscReal, PetscReal, PetscReal *, PetscReal *);
49: PetscErrorCode (*load)(TS, PetscViewer);
50: PetscErrorCode (*rollback)(TS);
51: PetscErrorCode (*getstages)(TS, PetscInt *, Vec *[]);
52: PetscErrorCode (*adjointstep)(TS);
53: PetscErrorCode (*adjointsetup)(TS);
54: PetscErrorCode (*adjointreset)(TS);
55: PetscErrorCode (*adjointintegral)(TS);
56: PetscErrorCode (*forwardsetup)(TS);
57: PetscErrorCode (*forwardreset)(TS);
58: PetscErrorCode (*forwardstep)(TS);
59: PetscErrorCode (*forwardintegral)(TS);
60: PetscErrorCode (*forwardgetstages)(TS, PetscInt *, Mat *[]);
61: PetscErrorCode (*getsolutioncomponents)(TS, PetscInt *, Vec *);
62: PetscErrorCode (*getauxsolution)(TS, Vec *);
63: PetscErrorCode (*gettimeerror)(TS, PetscInt, Vec *);
64: PetscErrorCode (*settimeerror)(TS, Vec);
65: PetscErrorCode (*startingmethod)(TS);
66: PetscErrorCode (*initcondition)(TS, Vec);
67: PetscErrorCode (*exacterror)(TS, Vec, Vec);
68: PetscErrorCode (*resizeregister)(TS, PetscBool);
69: };
71: /*
72: TSEvent - Abstract object to handle event monitoring
73: */
74: typedef struct _n_TSEvent *TSEvent;
76: typedef struct _TSTrajectoryOps *TSTrajectoryOps;
78: struct _TSTrajectoryOps {
79: PetscErrorCode (*view)(TSTrajectory, PetscViewer);
80: PetscErrorCode (*reset)(TSTrajectory);
81: PetscErrorCode (*destroy)(TSTrajectory);
82: PetscErrorCode (*set)(TSTrajectory, TS, PetscInt, PetscReal, Vec);
83: PetscErrorCode (*get)(TSTrajectory, TS, PetscInt, PetscReal *);
84: PetscErrorCode (*setfromoptions)(TSTrajectory, PetscOptionItems *);
85: PetscErrorCode (*setup)(TSTrajectory, TS);
86: };
88: /* TSHistory is an helper object that allows inquiring
89: the TSTrajectory by time and not by the step number only */
90: typedef struct _n_TSHistory *TSHistory;
92: struct _p_TSTrajectory {
93: PETSCHEADER(struct _TSTrajectoryOps);
94: TSHistory tsh; /* associates times to unique step ids */
95: /* stores necessary data to reconstruct states and derivatives via Lagrangian interpolation */
96: struct {
97: PetscInt order; /* interpolation order */
98: Vec *W; /* work vectors */
99: PetscScalar *L; /* workspace for Lagrange basis */
100: PetscReal *T; /* Lagrange times (stored) */
101: Vec *WW; /* just an array of pointers */
102: PetscBool *TT; /* workspace for Lagrange */
103: PetscReal *TW; /* Lagrange times (workspace) */
105: /* caching */
106: PetscBool caching;
107: struct {
108: PetscObjectId id;
109: PetscObjectState state;
110: PetscReal time;
111: PetscInt step;
112: } Ucached;
113: struct {
114: PetscObjectId id;
115: PetscObjectState state;
116: PetscReal time;
117: PetscInt step;
118: } Udotcached;
119: } lag;
120: Vec U, Udot; /* used by TSTrajectory{Get|Restore}UpdatedHistoryVecs */
121: PetscBool usehistory; /* whether to use TSHistory */
122: PetscBool solution_only; /* whether we dump just the solution or also the stages */
123: PetscBool adjoint_solve_mode; /* whether we will use the Trajectory inside a TSAdjointSolve() or not */
124: PetscViewer monitor;
125: PetscInt setupcalled; /* true if setup has been called */
126: PetscInt recomps; /* counter for recomputations in the adjoint run */
127: PetscInt diskreads, diskwrites; /* counters for disk checkpoint reads and writes */
128: char **names; /* the name of each variable; each process has only the local names */
129: PetscBool keepfiles; /* keep the files generated during the run after the run is complete */
130: char *dirname, *filetemplate; /* directory name and file name template for disk checkpoints */
131: char *dirfiletemplate; /* complete directory and file name template for disk checkpoints */
132: PetscErrorCode (*transform)(void *, Vec, Vec *);
133: PetscErrorCode (*transformdestroy)(void *);
134: void *transformctx;
135: void *data;
136: };
138: typedef struct _TS_RHSSplitLink *TS_RHSSplitLink;
139: struct _TS_RHSSplitLink {
140: TS ts;
141: char *splitname;
142: IS is;
143: TS_RHSSplitLink next;
144: PetscLogEvent event;
145: };
147: typedef struct _TS_EvaluationTimes *TSEvaluationTimes;
148: struct _TS_EvaluationTimes {
149: PetscInt num_time_points; /* number of time points */
150: PetscReal *time_points; /* array of the time span */
151: PetscReal reltol; /* relative tolerance for span point detection */
152: PetscReal abstol; /* absolute tolerance for span point detection */
153: PetscReal worktol; /* the ultimate tolerance (variable), maintained within a single TS time step for consistency */
154: PetscInt time_point_idx; /* index of the time_point to be reached next */
155: PetscInt sol_ctr; /* counter of the time points that have been reached */
156: Vec *sol_vecs; /* array of the solutions at the specified time points */
157: PetscReal *sol_times; /* array of times that sol_vecs was taken at */
158: };
160: struct _p_TS {
161: PETSCHEADER(struct _TSOps);
162: TSProblemType problem_type;
163: TSEquationType equation_type;
165: DM dm;
166: Vec vec_sol; /* solution vector in first and second order equations */
167: Vec vec_sol0; /* solution vector at the beginning of the step */
168: Vec vec_dot; /* time derivative vector in second order equations */
169: TSAdapt adapt;
170: TSAdaptType default_adapt_type;
171: TSEvent event;
173: /* ---------------- Resize ---------------------*/
174: PetscBool resizerollback;
175: PetscObjectList resizetransferobjs;
177: /* ---------------- User (or PETSc) Provided stuff ---------------------*/
178: PetscErrorCode (*monitor[MAXTSMONITORS])(TS, PetscInt, PetscReal, Vec, void *);
179: PetscCtxDestroyFn *monitordestroy[MAXTSMONITORS];
180: void *monitorcontext[MAXTSMONITORS];
181: PetscInt numbermonitors;
182: PetscErrorCode (*adjointmonitor[MAXTSMONITORS])(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *, void *);
183: PetscCtxDestroyFn *adjointmonitordestroy[MAXTSMONITORS];
184: void *adjointmonitorcontext[MAXTSMONITORS];
185: PetscInt numberadjointmonitors;
186: PetscInt monitorFrequency; /* Number of timesteps between monitor output */
188: PetscErrorCode (*prestep)(TS);
189: PetscErrorCode (*prestage)(TS, PetscReal);
190: PetscErrorCode (*poststage)(TS, PetscReal, PetscInt, Vec *);
191: PetscErrorCode (*postevaluate)(TS);
192: PetscErrorCode (*poststep)(TS);
193: PetscErrorCode (*functiondomainerror)(TS, PetscReal, Vec, PetscBool *);
194: PetscErrorCode (*resizesetup)(TS, PetscInt, PetscReal, Vec, PetscBool *, void *);
195: PetscErrorCode (*resizetransfer)(TS, PetscInt, Vec[], Vec[], void *);
196: void *resizectx;
198: /* ---------------------- Sensitivity Analysis support -----------------*/
199: TSTrajectory trajectory; /* All solutions are kept here for the entire time integration process */
200: Vec *vecs_sensi; /* one vector for each cost function */
201: Vec *vecs_sensip;
202: PetscInt numcost; /* number of cost functions */
203: Vec vec_costintegral;
204: PetscInt adjointsetupcalled;
205: PetscInt adjoint_steps;
206: PetscInt adjoint_max_steps;
207: PetscBool adjoint_solve; /* immediately call TSAdjointSolve() after TSSolve() is complete */
208: PetscBool costintegralfwd; /* cost integral is evaluated in the forward run if true */
209: Vec vec_costintegrand; /* workspace for Adjoint computations */
210: Mat Jacp, Jacprhs;
211: void *ijacobianpctx, *rhsjacobianpctx;
212: void *costintegrandctx;
213: Vec *vecs_drdu;
214: Vec *vecs_drdp;
215: Vec vec_drdu_col, vec_drdp_col;
217: /* first-order adjoint */
218: PetscErrorCode (*rhsjacobianp)(TS, PetscReal, Vec, Mat, void *);
219: PetscErrorCode (*ijacobianp)(TS, PetscReal, Vec, Vec, PetscReal, Mat, void *);
220: PetscErrorCode (*costintegrand)(TS, PetscReal, Vec, Vec, void *);
221: PetscErrorCode (*drdufunction)(TS, PetscReal, Vec, Vec *, void *);
222: PetscErrorCode (*drdpfunction)(TS, PetscReal, Vec, Vec *, void *);
224: /* second-order adjoint */
225: Vec *vecs_sensi2;
226: Vec *vecs_sensi2p;
227: Vec vec_dir; /* directional vector for optimization */
228: Vec *vecs_fuu, *vecs_guu;
229: Vec *vecs_fup, *vecs_gup;
230: Vec *vecs_fpu, *vecs_gpu;
231: Vec *vecs_fpp, *vecs_gpp;
232: void *ihessianproductctx, *rhshessianproductctx;
233: PetscErrorCode (*ihessianproduct_fuu)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *);
234: PetscErrorCode (*ihessianproduct_fup)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *);
235: PetscErrorCode (*ihessianproduct_fpu)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *);
236: PetscErrorCode (*ihessianproduct_fpp)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *);
237: PetscErrorCode (*rhshessianproduct_guu)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *);
238: PetscErrorCode (*rhshessianproduct_gup)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *);
239: PetscErrorCode (*rhshessianproduct_gpu)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *);
240: PetscErrorCode (*rhshessianproduct_gpp)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *);
242: /* specific to forward sensitivity analysis */
243: Mat mat_sensip; /* matrix storing forward sensitivities */
244: Vec vec_sensip_col; /* space for a column of the sensip matrix */
245: Vec *vecs_integral_sensip; /* one vector for each integral */
246: PetscInt num_parameters;
247: PetscInt num_initialvalues;
248: void *vecsrhsjacobianpctx;
249: PetscInt forwardsetupcalled;
250: PetscBool forward_solve;
251: PetscErrorCode (*vecsrhsjacobianp)(TS, PetscReal, Vec, Vec *, void *);
253: /* ---------------------- IMEX support ---------------------------------*/
254: /* These extra slots are only used when the user provides both Implicit and RHS */
255: Mat Arhs; /* Right hand side matrix */
256: Mat Brhs; /* Right hand side preconditioning matrix */
257: Vec Frhs; /* Right hand side function value */
259: /* This is a general caching scheme to avoid recomputing the Jacobian at a place that has been previously been evaluated.
260: * The present use case is that TSComputeRHSFunctionLinear() evaluates the Jacobian once and we don't want it to be immeditely re-evaluated.
261: */
262: struct {
263: PetscReal time; /* The time at which the matrices were last evaluated */
264: PetscObjectId Xid; /* Unique ID of solution vector at which the Jacobian was last evaluated */
265: PetscObjectState Xstate; /* State of the solution vector */
266: MatStructure mstructure; /* The structure returned */
267: /* Flag to unshift Jacobian before calling the IJacobian or RHSJacobian functions. This is useful
268: * if the user would like to reuse (part of) the Jacobian from the last evaluation. */
269: PetscBool reuse;
270: PetscReal scale, shift;
271: } rhsjacobian;
273: struct {
274: PetscReal shift; /* The derivative of the lhs wrt to Xdot */
275: } ijacobian;
277: MatStructure axpy_pattern; /* information about the nonzero pattern of the RHS Jacobian in reference to the implicit Jacobian */
278: /* --------------------Nonlinear Iteration------------------------------*/
279: SNES snes;
280: PetscBool usessnes; /* Flag set by each TSType to indicate if the type actually uses a SNES;
281: this works around the design flaw that a SNES is ALWAYS created with TS even when it is not needed.*/
282: PetscInt ksp_its; /* total number of linear solver iterations */
283: PetscInt snes_its; /* total number of nonlinear solver iterations */
284: PetscInt num_snes_failures;
285: PetscInt max_snes_failures;
287: /* --- Logging --- */
288: PetscInt ifuncs, rhsfuncs, ijacs, rhsjacs;
290: /* --- Data that is unique to each particular solver --- */
291: PetscInt setupcalled; /* true if setup has been called */
292: void *data; /* implementationspecific data */
293: void *ctx; /* user context */
295: PetscBool steprollback; /* flag to indicate that the step was rolled back */
296: PetscBool steprestart; /* flag to indicate that the timestepper has to discard any history and restart */
297: PetscBool stepresize; /* flag to indicate that the discretization was resized */
298: PetscInt steps; /* steps taken so far in all successive calls to TSSolve() */
299: PetscReal ptime; /* time at the start of the current step (stage time is internal if it exists) */
300: PetscReal time_step; /* current time increment */
301: PetscReal time_step0; /* proposed time increment at the beginning of the step */
302: PetscReal ptime_prev; /* time at the start of the previous step */
303: PetscReal ptime_prev_rollback; /* time at the start of the 2nd previous step to recover from rollback */
304: PetscReal solvetime; /* time at the conclusion of TSSolve() */
305: PetscBool stifflyaccurate; /* flag to indicate that the method is stiffly accurate */
307: TSConvergedReason reason;
308: PetscBool errorifstepfailed;
309: PetscInt reject, max_reject;
310: TSExactFinalTimeOption exact_final_time;
312: PetscObjectParameterDeclare(PetscReal, rtol); /* Relative and absolute tolerance for local truncation error */
313: PetscObjectParameterDeclare(PetscReal, atol);
314: PetscObjectParameterDeclare(PetscReal, max_time); /* max time allowed */
315: PetscObjectParameterDeclare(PetscInt, max_steps); /* max number of steps */
316: Vec vatol, vrtol; /* Relative and absolute tolerance in vector form */
317: PetscReal cfltime, cfltime_local;
319: PetscBool testjacobian;
320: PetscBool testjacobiantranspose;
321: /* ------------------- Default work-area management ------------------ */
322: PetscInt nwork;
323: Vec *work;
325: /* ---------------------- RHS splitting support ---------------------------------*/
326: PetscInt num_rhs_splits;
327: TS_RHSSplitLink tsrhssplit;
328: PetscBool use_splitrhsfunction;
329: SNES snesrhssplit;
331: /* ---------------------- Quadrature integration support ---------------------------------*/
332: TS quadraturets;
334: /* ---------------------- Time span support ---------------------------------*/
335: TSEvaluationTimes eval_times;
336: };
338: struct _TSAdaptOps {
339: PetscErrorCode (*choose)(TSAdapt, TS, PetscReal, PetscInt *, PetscReal *, PetscBool *, PetscReal *, PetscReal *, PetscReal *);
340: PetscErrorCode (*destroy)(TSAdapt);
341: PetscErrorCode (*reset)(TSAdapt);
342: PetscErrorCode (*view)(TSAdapt, PetscViewer);
343: PetscErrorCode (*setfromoptions)(TSAdapt, PetscOptionItems *);
344: PetscErrorCode (*load)(TSAdapt, PetscViewer);
345: };
347: struct _p_TSAdapt {
348: PETSCHEADER(struct _TSAdaptOps);
349: void *data;
350: PetscErrorCode (*checkstage)(TSAdapt, TS, PetscReal, Vec, PetscBool *);
351: struct {
352: PetscInt n; /* number of candidate schemes, including the one currently in use */
353: PetscBool inuse_set; /* the current scheme has been set */
354: const char *name[16]; /* name of the scheme */
355: PetscInt order[16]; /* classical order of each scheme */
356: PetscInt stageorder[16]; /* stage order of each scheme */
357: PetscReal ccfl[16]; /* stability limit relative to explicit Euler */
358: PetscReal cost[16]; /* relative measure of the amount of work required for each scheme */
359: } candidates;
360: PetscBool always_accept;
361: PetscReal safety; /* safety factor relative to target error/stability goal */
362: PetscReal reject_safety; /* extra safety factor if the last step was rejected */
363: PetscReal clip[2]; /* admissible time step decrease/increase factors */
364: PetscReal dt_min, dt_max; /* admissible minimum and maximum time step */
365: PetscReal ignore_max; /* minimum value of the solution to be considered by the adaptor */
366: PetscBool glee_use_local; /* GLEE adaptor uses global or local error */
367: PetscReal scale_solve_failed; /* scale step by this factor if solver (linear or nonlinear) fails. */
368: PetscReal matchstepfac[2]; /* factors to control the behaviour of matchstep */
369: NormType wnormtype;
370: PetscViewer monitor;
371: PetscInt timestepjustdecreased_delay; /* number of timesteps after a decrease in the timestep before the timestep can be increased */
372: PetscInt timestepjustdecreased;
373: PetscReal dt_eval_times_cached; /* time step before hitting a TS evaluation time point */
374: };
376: typedef struct _p_DMTS *DMTS;
377: typedef struct _DMTSOps *DMTSOps;
378: struct _DMTSOps {
379: TSRHSFunctionFn *rhsfunction;
380: TSRHSJacobianFn *rhsjacobian;
382: TSIFunctionFn *ifunction;
383: PetscErrorCode (*ifunctionview)(void *, PetscViewer);
384: PetscErrorCode (*ifunctionload)(void **, PetscViewer);
386: TSIJacobianFn *ijacobian;
387: PetscErrorCode (*ijacobianview)(void *, PetscViewer);
388: PetscErrorCode (*ijacobianload)(void **, PetscViewer);
390: TSI2FunctionFn *i2function;
391: TSI2JacobianFn *i2jacobian;
393: TSTransientVariableFn *transientvar;
395: TSSolutionFn *solution;
396: TSForcingFn *forcing;
398: PetscErrorCode (*destroy)(DMTS);
399: PetscErrorCode (*duplicate)(DMTS, DMTS);
400: };
402: /*S
403: DMTS - Object held by a `DM` that contains all the callback functions and their contexts needed by a `TS`
405: Level: developer
407: Notes:
408: Users provide callback functions and their contexts to `TS` using, for example, `TSSetIFunction()`. These values are stored
409: in a `DMTS` that is contained in the `DM` associated with the `TS`. If no `DM` was provided by
410: the user with `TSSetDM()` it is automatically created by `TSGetDM()` with `DMShellCreate()`.
412: Users very rarely need to worked directly with the `DMTS` object, rather they work with the `TS` and the `DM` they created
414: Multiple `DM` can share a single `DMTS`, often each `DM` is associated with
415: a grid refinement level. `DMGetDMTS()` returns the `DMTS` associated with a `DM`. `DMGetDMTSWrite()` returns a unique
416: `DMTS` that is only associated with the current `DM`, making a copy of the shared `DMTS` if needed (copy-on-write).
418: See `DMKSP` for details on why there is a needed for `DMTS` instead of simply storing the user callbacks directly in the `DM` or the `TS`
420: Developer Note:
421: The original `dm` inside the `DMTS` is NOT reference counted (to prevent a reference count loop between a `DM` and a `DMSNES`).
422: The `DM` on which this context was first created is cached here to implement one-way
423: copy-on-write. When `DMGetDMTSWrite()` sees a request using a different `DM`, it makes a copy of the `DMTS`.
425: .seealso: [](ch_ts), `TSCreate()`, `DM`, `DMGetDMTSWrite()`, `DMGetDMTS()`, `TSSetIFunction()`, `DMTSSetRHSFunctionContextDestroy()`,
426: `DMTSSetRHSJacobian()`, `DMTSGetRHSJacobian()`, `DMTSSetRHSJacobianContextDestroy()`, `DMTSSetIFunction()`, `DMTSGetIFunction()`,
427: `DMTSSetIFunctionContextDestroy()`, `DMTSSetIJacobian()`, `DMTSGetIJacobian()`, `DMTSSetIJacobianContextDestroy()`,
428: `DMTSSetI2Function()`, `DMTSGetI2Function()`, `DMTSSetI2FunctionContextDestroy()`, `DMTSSetI2Jacobian()`,
429: `DMTSGetI2Jacobian()`, `DMTSSetI2JacobianContextDestroy()`, `DMKSP`, `DMSNES`
430: S*/
431: struct _p_DMTS {
432: PETSCHEADER(struct _DMTSOps);
433: PetscContainer rhsfunctionctxcontainer;
434: PetscContainer rhsjacobianctxcontainer;
436: PetscContainer ifunctionctxcontainer;
437: PetscContainer ijacobianctxcontainer;
439: PetscContainer i2functionctxcontainer;
440: PetscContainer i2jacobianctxcontainer;
442: void *transientvarctx;
444: void *solutionctx;
445: void *forcingctx;
447: void *data;
449: /* See the developer note for DMTS above */
450: DM originaldm;
451: };
453: PETSC_EXTERN PetscErrorCode DMTSUnsetRHSFunctionContext_Internal(DM);
454: PETSC_EXTERN PetscErrorCode DMTSUnsetRHSJacobianContext_Internal(DM);
455: PETSC_EXTERN PetscErrorCode DMTSUnsetIFunctionContext_Internal(DM);
456: PETSC_EXTERN PetscErrorCode DMTSUnsetIJacobianContext_Internal(DM);
457: PETSC_EXTERN PetscErrorCode DMTSUnsetI2FunctionContext_Internal(DM);
458: PETSC_EXTERN PetscErrorCode DMTSUnsetI2JacobianContext_Internal(DM);
460: PETSC_EXTERN PetscErrorCode DMGetDMTS(DM, DMTS *);
461: PETSC_EXTERN PetscErrorCode DMGetDMTSWrite(DM, DMTS *);
462: PETSC_EXTERN PetscErrorCode DMCopyDMTS(DM, DM);
463: PETSC_EXTERN PetscErrorCode DMTSView(DMTS, PetscViewer);
464: PETSC_EXTERN PetscErrorCode DMTSLoad(DMTS, PetscViewer);
465: PETSC_EXTERN PetscErrorCode DMTSCopy(DMTS, DMTS);
467: struct _n_TSEvent {
468: PetscReal *fvalue_prev; /* value of indicator function at the left end-point of the event interval */
469: PetscReal *fvalue; /* value of indicator function at the current point */
470: PetscReal *fvalue_right; /* value of indicator function at the right end-point of the event interval */
471: PetscInt *fsign_prev; /* sign of indicator function at the left end-point of the event interval */
472: PetscInt *fsign; /* sign of indicator function at the current point */
473: PetscInt *fsign_right; /* sign of indicator function at the right end-point of the event interval */
474: PetscReal ptime_prev; /* time at the previous point (left end-point of the event interval) */
475: PetscReal ptime_right; /* time at the right end-point of the event interval */
476: PetscReal ptime_cache; /* the point visited by the TS before the event interval was detected; cached - to reuse if necessary */
477: PetscReal timestep_cache; /* time step considered by the TS before the event interval was detected; cached - to reuse if necessary */
478: PetscInt *side; /* upon bracket subdivision, indicates which sub-bracket is taken further, -1 -> left one, +1 -> right one, +2 -> neither, 0 -> zero-crossing located */
479: PetscInt *side_prev; /* counts the repeating previous side's (with values: -n <=> '-1'*n; +n <=> '+1'*n); used in the Anderson-Bjorck iteration */
480: PetscReal timestep_postevent; /* first time step after the event; can be PETSC_DECIDE */
481: PetscReal timestep_2nd_postevent; /* second time step after the event; can be PETSC_DECIDE */
482: PetscReal timestep_min; /* minimum time step */
483: PetscBool *justrefined_AB; /* this flag shows if the given indicator function i = [0..nevents) participated in Anderson-Bjorck process in the last iteration of TSEventHandler() */
484: PetscReal *gamma_AB; /* cumulative scaling factor for the Anderson-Bjorck iteration */
485: PetscErrorCode (*indicator)(TS, PetscReal, Vec, PetscReal *, void *); /* this callback defines the user function(s) whose sign changes indicate events */
486: PetscErrorCode (*postevent)(TS, PetscInt, PetscInt[], PetscReal, Vec, PetscBool, void *); /* user post-event callback */
487: void *ctx; /* user context for indicator and postevent callbacks */
488: PetscInt *direction; /* zero crossing direction to trigger the event: +1 -> going positive, -1 -> going negative, 0 -> any */
489: PetscBool *terminate; /* 1 -> terminate time stepping on event location, 0 -> continue */
490: PetscInt nevents; /* number of events (indicator functions) to handle on the current MPI process */
491: PetscInt nevents_zero; /* number of events triggered */
492: PetscInt *events_zero; /* list of the events triggered */
493: PetscReal *vtol; /* array of tolerances for the indicator function zero check */
494: PetscInt iterctr; /* iteration counter: used both for reporting and as a status indicator */
495: PetscBool processing; /* this flag shows if the event-resolving iterations are in progress, or the post-event dt handling is in progress */
496: PetscBool revisit_right; /* [sync] "revisit the bracket's right end", if true, then fvalue(s) are not calculated, but are taken from fvalue_right(s) */
497: PetscViewer monitor;
498: /* Struct to record the events */
499: struct {
500: PetscInt ctr; /* Recorder counter */
501: PetscReal *time; /* Event times */
502: PetscInt *stepnum; /* Step numbers */
503: PetscInt *nevents; /* Number of events occurring at the event times */
504: PetscInt **eventidx; /* Local indices of the events in the event list */
505: } recorder;
506: PetscInt recsize; /* Size of recorder stack */
507: PetscInt refct; /* Reference count */
508: };
510: PETSC_EXTERN PetscErrorCode TSEventInitialize(TSEvent, TS, PetscReal, Vec);
511: PETSC_EXTERN PetscErrorCode TSEventDestroy(TSEvent *);
512: PETSC_EXTERN PetscErrorCode TSEventHandler(TS);
513: PETSC_EXTERN PetscErrorCode TSAdjointEventHandler(TS);
515: PETSC_EXTERN PetscLogEvent TS_AdjointStep;
516: PETSC_EXTERN PetscLogEvent TS_Step;
517: PETSC_EXTERN PetscLogEvent TS_PseudoComputeTimeStep;
518: PETSC_EXTERN PetscLogEvent TS_FunctionEval;
519: PETSC_EXTERN PetscLogEvent TS_JacobianEval;
520: PETSC_EXTERN PetscLogEvent TS_ForwardStep;
522: typedef enum {
523: TS_STEP_INCOMPLETE, /* vec_sol, ptime, etc point to beginning of step */
524: TS_STEP_PENDING, /* vec_sol advanced, but step has not been accepted yet */
525: TS_STEP_COMPLETE /* step accepted and ptime, steps, etc have been advanced */
526: } TSStepStatus;
528: struct _n_TSMonitorLGCtx {
529: PetscDrawLG lg;
530: PetscBool semilogy;
531: PetscInt howoften; /* when > 0 uses step % howoften, when negative only final solution plotted */
532: PetscInt ksp_its, snes_its;
533: char **names;
534: char **displaynames;
535: PetscInt ndisplayvariables;
536: PetscInt *displayvariables;
537: PetscReal *displayvalues;
538: PetscErrorCode (*transform)(void *, Vec, Vec *);
539: PetscCtxDestroyFn *transformdestroy;
540: void *transformctx;
541: };
543: struct _n_TSMonitorSPCtx {
544: PetscDrawSP sp;
545: PetscInt howoften; /* when > 0 uses step % howoften, when negative only final solution plotted */
546: PetscInt retain; /* Retain n points plotted to show trajectories, or -1 for all points */
547: PetscBool phase; /* Plot in phase space rather than coordinate space */
548: PetscBool multispecies; /* Change scatter point color based on species */
549: PetscInt ksp_its, snes_its;
550: };
552: struct _n_TSMonitorHGCtx {
553: PetscDrawHG *hg;
554: PetscInt howoften; /* when > 0 uses step % howoften, when negative only final solution plotted */
555: PetscInt Ns; /* The number of species to histogram */
556: PetscBool velocity; /* Plot in velocity space rather than coordinate space */
557: };
559: struct _n_TSMonitorEnvelopeCtx {
560: Vec max, min;
561: };
563: /*
564: Checks if the user provide a TSSetIFunction() but an explicit method is called; generate an error in that case
565: */
566: static inline PetscErrorCode TSCheckImplicitTerm(TS ts)
567: {
568: TSIFunctionFn *ifunction;
569: DM dm;
571: PetscFunctionBegin;
572: PetscCall(TSGetDM(ts, &dm));
573: PetscCall(DMTSGetIFunction(dm, &ifunction, NULL));
574: PetscCheck(!ifunction, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_INCOMP, "You are attempting to use an explicit ODE integrator but provided an implicit function definition with TSSetIFunction()");
575: PetscFunctionReturn(PETSC_SUCCESS);
576: }
578: PETSC_EXTERN PetscErrorCode TSGetRHSMats_Private(TS, Mat *, Mat *);
579: /* this is declared here as TSHistory is not public */
580: PETSC_EXTERN PetscErrorCode TSAdaptHistorySetTSHistory(TSAdapt, TSHistory, PetscBool);
582: PETSC_INTERN PetscErrorCode TSTrajectoryReconstruct_Private(TSTrajectory, TS, PetscReal, Vec, Vec);
583: PETSC_INTERN PetscErrorCode TSTrajectorySetUp_Basic(TSTrajectory, TS);
585: PETSC_EXTERN PetscLogEvent TSTrajectory_Set;
586: PETSC_EXTERN PetscLogEvent TSTrajectory_Get;
587: PETSC_EXTERN PetscLogEvent TSTrajectory_GetVecs;
588: PETSC_EXTERN PetscLogEvent TSTrajectory_SetUp;
589: PETSC_EXTERN PetscLogEvent TSTrajectory_DiskWrite;
590: PETSC_EXTERN PetscLogEvent TSTrajectory_DiskRead;
592: struct _n_TSMonitorDrawCtx {
593: PetscViewer viewer;
594: Vec initialsolution;
595: PetscBool showinitial;
596: PetscInt howoften; /* when > 0 uses step % howoften, when negative only final solution plotted */
597: PetscBool showtimestepandtime;
598: };
600: struct _n_TSMonitorVTKCtx {
601: char *filenametemplate;
602: PetscInt interval; /* when > 0 uses step % interval, when negative only final solution plotted */
603: };