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