Actual source code: tsimpl.h
petsc-3.12.1 2019-10-22
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);
31: typedef struct _TSOps *TSOps;
33: struct _TSOps {
34: PetscErrorCode (*snesfunction)(SNES,Vec,Vec,TS);
35: PetscErrorCode (*snesjacobian)(SNES,Vec,Mat,Mat,TS);
36: PetscErrorCode (*setup)(TS);
37: PetscErrorCode (*step)(TS);
38: PetscErrorCode (*solve)(TS);
39: PetscErrorCode (*interpolate)(TS,PetscReal,Vec);
40: PetscErrorCode (*evaluatewlte)(TS,NormType,PetscInt*,PetscReal*);
41: PetscErrorCode (*evaluatestep)(TS,PetscInt,Vec,PetscBool*);
42: PetscErrorCode (*setfromoptions)(PetscOptionItems*,TS);
43: PetscErrorCode (*destroy)(TS);
44: PetscErrorCode (*view)(TS,PetscViewer);
45: PetscErrorCode (*reset)(TS);
46: PetscErrorCode (*linearstability)(TS,PetscReal,PetscReal,PetscReal*,PetscReal*);
47: PetscErrorCode (*load)(TS,PetscViewer);
48: PetscErrorCode (*rollback)(TS);
49: PetscErrorCode (*getstages)(TS,PetscInt*,Vec**);
50: PetscErrorCode (*adjointstep)(TS);
51: PetscErrorCode (*adjointsetup)(TS);
52: PetscErrorCode (*adjointreset)(TS);
53: PetscErrorCode (*adjointintegral)(TS);
54: PetscErrorCode (*forwardsetup)(TS);
55: PetscErrorCode (*forwardreset)(TS);
56: PetscErrorCode (*forwardstep)(TS);
57: PetscErrorCode (*forwardintegral)(TS);
58: PetscErrorCode (*forwardgetstages)(TS,PetscInt*,Mat**);
59: PetscErrorCode (*getsolutioncomponents)(TS,PetscInt*,Vec*);
60: PetscErrorCode (*getauxsolution)(TS,Vec*);
61: PetscErrorCode (*gettimeerror)(TS,PetscInt,Vec*);
62: PetscErrorCode (*settimeerror)(TS,Vec);
63: PetscErrorCode (*startingmethod) (TS);
64: };
66: /*
67: TSEvent - Abstract object to handle event monitoring
68: */
69: typedef struct _n_TSEvent *TSEvent;
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)(PetscOptionItems*,TSTrajectory);
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: PetscInt 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: PetscErrorCode (*transformdestroy)(void*);
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: struct _p_TS {
143: PETSCHEADER(struct _TSOps);
144: TSProblemType problem_type;
145: TSEquationType equation_type;
147: DM dm;
148: Vec vec_sol; /* solution vector in first and second order equations */
149: Vec vec_dot; /* time derivative vector in second order equations */
150: TSAdapt adapt;
151: TSAdaptType default_adapt_type;
152: TSEvent event;
154: /* ---------------- User (or PETSc) Provided stuff ---------------------*/
155: PetscErrorCode (*monitor[MAXTSMONITORS])(TS,PetscInt,PetscReal,Vec,void*);
156: PetscErrorCode (*monitordestroy[MAXTSMONITORS])(void**);
157: void *monitorcontext[MAXTSMONITORS];
158: PetscInt numbermonitors;
159: PetscErrorCode (*adjointmonitor[MAXTSMONITORS])(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*);
160: PetscErrorCode (*adjointmonitordestroy[MAXTSMONITORS])(void**);
161: void *adjointmonitorcontext[MAXTSMONITORS];
162: PetscInt numberadjointmonitors;
164: PetscErrorCode (*prestep)(TS);
165: PetscErrorCode (*prestage)(TS,PetscReal);
166: PetscErrorCode (*poststage)(TS,PetscReal,PetscInt,Vec*);
167: PetscErrorCode (*postevaluate)(TS);
168: PetscErrorCode (*poststep)(TS);
169: PetscErrorCode (*functiondomainerror)(TS,PetscReal,Vec,PetscBool*);
171: /* ---------------------- Sensitivity Analysis support -----------------*/
172: TSTrajectory trajectory; /* All solutions are kept here for the entire time integration process */
173: Vec *vecs_sensi; /* one vector for each cost function */
174: Vec *vecs_sensip;
175: PetscInt numcost; /* number of cost functions */
176: Vec vec_costintegral;
177: PetscInt adjointsetupcalled;
178: PetscInt adjoint_steps;
179: PetscInt adjoint_max_steps;
180: PetscBool adjoint_solve; /* immediately call TSAdjointSolve() after TSSolve() is complete */
181: PetscBool costintegralfwd; /* cost integral is evaluated in the forward run if true */
182: Vec vec_costintegrand; /* workspace for Adjoint computations */
183: Mat Jacp,Jacprhs;
184: void *ijacobianpctx,*rhsjacobianpctx;
185: void *costintegrandctx;
186: Vec *vecs_drdu;
187: Vec *vecs_drdp;
188: Vec vec_drdu_col,vec_drdp_col;
190: /* first-order adjoint */
191: PetscErrorCode (*rhsjacobianp)(TS,PetscReal,Vec,Mat,void*);
192: PetscErrorCode (*ijacobianp)(TS,PetscReal,Vec,Vec,PetscReal,Mat,void*);
193: PetscErrorCode (*costintegrand)(TS,PetscReal,Vec,Vec,void*);
194: PetscErrorCode (*drdufunction)(TS,PetscReal,Vec,Vec*,void*);
195: PetscErrorCode (*drdpfunction)(TS,PetscReal,Vec,Vec*,void*);
197: /* second-order adjoint */
198: Vec *vecs_sensi2;
199: Vec *vecs_sensi2p;
200: Vec vec_dir; /* directional vector for optimization */
201: Vec *vecs_fuu,*vecs_guu;
202: Vec *vecs_fup,*vecs_gup;
203: Vec *vecs_fpu,*vecs_gpu;
204: Vec *vecs_fpp,*vecs_gpp;
205: void *ihessianproductctx,*rhshessianproductctx;
206: PetscErrorCode (*ihessianproduct_fuu)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*);
207: PetscErrorCode (*ihessianproduct_fup)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*);
208: PetscErrorCode (*ihessianproduct_fpu)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*);
209: PetscErrorCode (*ihessianproduct_fpp)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*);
210: PetscErrorCode (*rhshessianproduct_guu)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*);
211: PetscErrorCode (*rhshessianproduct_gup)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*);
212: PetscErrorCode (*rhshessianproduct_gpu)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*);
213: PetscErrorCode (*rhshessianproduct_gpp)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*);
215: /* specific to forward sensitivity analysis */
216: Mat mat_sensip; /* matrix storing forward sensitivities */
217: Vec vec_sensip_col; /* space for a column of the sensip matrix */
218: Vec *vecs_integral_sensip; /* one vector for each integral */
219: PetscInt num_parameters;
220: PetscInt num_initialvalues;
221: void *vecsrhsjacobianpctx;
222: PetscInt forwardsetupcalled;
223: PetscBool forward_solve;
224: PetscErrorCode (*vecsrhsjacobianp)(TS,PetscReal,Vec,Vec*,void*);
226: /* ---------------------- IMEX support ---------------------------------*/
227: /* These extra slots are only used when the user provides both Implicit and RHS */
228: Mat Arhs; /* Right hand side matrix */
229: Mat Brhs; /* Right hand side preconditioning matrix */
230: Vec Frhs; /* Right hand side function value */
232: /* This is a general caching scheme to avoid recomputing the Jacobian at a place that has been previously been evaluated.
233: * The present use case is that TSComputeRHSFunctionLinear() evaluates the Jacobian once and we don't want it to be immeditely re-evaluated.
234: */
235: struct {
236: PetscReal time; /* The time at which the matrices were last evaluated */
237: PetscObjectId Xid; /* Unique ID of solution vector at which the Jacobian was last evaluated */
238: PetscObjectState Xstate; /* State of the solution vector */
239: MatStructure mstructure; /* The structure returned */
240: /* Flag to unshift Jacobian before calling the IJacobian or RHSJacobian functions. This is useful
241: * if the user would like to reuse (part of) the Jacobian from the last evaluation. */
242: PetscBool reuse;
243: PetscReal scale,shift;
244: } rhsjacobian;
246: struct {
247: PetscReal shift; /* The derivative of the lhs wrt to Xdot */
248: } ijacobian;
250: /* --------------------Nonlinear Iteration------------------------------*/
251: SNES snes;
252: PetscBool usessnes; /* Flag set by each TSType to indicate if the type actually uses a SNES;
253: this works around the design flaw that a SNES is ALWAYS created with TS even when it is not needed.*/
254: PetscInt ksp_its; /* total number of linear solver iterations */
255: PetscInt snes_its; /* total number of nonlinear solver iterations */
256: PetscInt num_snes_failures;
257: PetscInt max_snes_failures;
259: /* --- Data that is unique to each particular solver --- */
260: PetscInt setupcalled; /* true if setup has been called */
261: void *data; /* implementationspecific data */
262: void *user; /* user context */
264: /* ------------------ Parameters -------------------------------------- */
265: PetscInt max_steps; /* max number of steps */
266: PetscReal max_time; /* max time allowed */
268: /* --------------------------------------------------------------------- */
270: PetscBool steprollback; /* flag to indicate that the step was rolled back */
271: PetscBool steprestart; /* flag to indicate that the timestepper has to discard any history and restart */
272: PetscInt steps; /* steps taken so far in all successive calls to TSSolve() */
273: PetscReal ptime; /* time at the start of the current step (stage time is internal if it exists) */
274: PetscReal time_step; /* current time increment */
275: PetscReal ptime_prev; /* time at the start of the previous step */
276: PetscReal ptime_prev_rollback; /* time at the start of the 2nd previous step to recover from rollback */
277: PetscReal solvetime; /* time at the conclusion of TSSolve() */
279: TSConvergedReason reason;
280: PetscBool errorifstepfailed;
281: PetscInt reject,max_reject;
282: TSExactFinalTimeOption exact_final_time;
284: PetscReal atol,rtol; /* Relative and absolute tolerance for local truncation error */
285: Vec vatol,vrtol; /* Relative and absolute tolerance in vector form */
286: PetscReal cfltime,cfltime_local;
288: PetscBool testjacobian;
289: PetscBool testjacobiantranspose;
290: /* ------------------- Default work-area management ------------------ */
291: PetscInt nwork;
292: Vec *work;
294: /* ---------------------- RHS splitting support ---------------------------------*/
295: PetscInt num_rhs_splits;
296: TS_RHSSplitLink tsrhssplit;
297: PetscBool use_splitrhsfunction;
299: /* ---------------------- Quadrature integration support ---------------------------------*/
300: TS quadraturets;
301: };
303: struct _TSAdaptOps {
304: PetscErrorCode (*choose)(TSAdapt,TS,PetscReal,PetscInt*,PetscReal*,PetscBool*,PetscReal*,PetscReal*,PetscReal*);
305: PetscErrorCode (*destroy)(TSAdapt);
306: PetscErrorCode (*reset)(TSAdapt);
307: PetscErrorCode (*view)(TSAdapt,PetscViewer);
308: PetscErrorCode (*setfromoptions)(PetscOptionItems*,TSAdapt);
309: PetscErrorCode (*load)(TSAdapt,PetscViewer);
310: };
312: struct _p_TSAdapt {
313: PETSCHEADER(struct _TSAdaptOps);
314: void *data;
315: PetscErrorCode (*checkstage)(TSAdapt,TS,PetscReal,Vec,PetscBool*);
316: struct {
317: PetscInt n; /* number of candidate schemes, including the one currently in use */
318: PetscBool inuse_set; /* the current scheme has been set */
319: const char *name[16]; /* name of the scheme */
320: PetscInt order[16]; /* classical order of each scheme */
321: PetscInt stageorder[16]; /* stage order of each scheme */
322: PetscReal ccfl[16]; /* stability limit relative to explicit Euler */
323: PetscReal cost[16]; /* relative measure of the amount of work required for each scheme */
324: } candidates;
325: PetscBool always_accept;
326: PetscReal safety; /* safety factor relative to target error/stability goal */
327: PetscReal reject_safety; /* extra safety factor if the last step was rejected */
328: PetscReal clip[2]; /* admissible time step decrease/increase factors */
329: PetscReal dt_min,dt_max; /* admissible minimum and maximum time step */
330: PetscReal ignore_max; /* minimum value of the solution to be considered by the adaptor */
331: PetscBool glee_use_local; /* GLEE adaptor uses global or local error */
332: PetscReal scale_solve_failed; /* scale step by this factor if solver (linear or nonlinear) fails. */
333: PetscReal matchstepfac[2]; /* factors to control the behaviour of matchstep */
334: NormType wnormtype;
335: PetscViewer monitor;
336: PetscInt timestepjustdecreased_delay; /* number of timesteps after a decrease in the timestep before the timestep can be increased */
337: PetscInt timestepjustdecreased;
338: };
340: typedef struct _p_DMTS *DMTS;
341: typedef struct _DMTSOps *DMTSOps;
342: struct _DMTSOps {
343: TSRHSFunction rhsfunction;
344: TSRHSJacobian rhsjacobian;
346: TSIFunction ifunction;
347: PetscErrorCode (*ifunctionview)(void*,PetscViewer);
348: PetscErrorCode (*ifunctionload)(void**,PetscViewer);
350: TSIJacobian ijacobian;
351: PetscErrorCode (*ijacobianview)(void*,PetscViewer);
352: PetscErrorCode (*ijacobianload)(void**,PetscViewer);
354: TSI2Function i2function;
355: TSI2Jacobian i2jacobian;
357: TSSolutionFunction solution;
358: TSForcingFunction forcing;
360: PetscErrorCode (*destroy)(DMTS);
361: PetscErrorCode (*duplicate)(DMTS,DMTS);
362: };
364: struct _p_DMTS {
365: PETSCHEADER(struct _DMTSOps);
366: void *rhsfunctionctx;
367: void *rhsjacobianctx;
369: void *ifunctionctx;
370: void *ijacobianctx;
372: void *i2functionctx;
373: void *i2jacobianctx;
375: void *solutionctx;
376: void *forcingctx;
378: void *data;
380: /* This is NOT reference counted. The DM on which this context was first created is cached here to implement one-way
381: * copy-on-write. When DMGetDMTSWrite() sees a request using a different DM, it makes a copy. Thus, if a user
382: * only interacts directly with one level, e.g., using TSSetIFunction(), then coarse levels of a multilevel item
383: * integrator are built, then the user changes the routine with another call to TSSetIFunction(), it automatically
384: * propagates to all the levels. If instead, they get out a specific level and set the function on that level,
385: * subsequent changes to the original level will no longer propagate to that level.
386: */
387: DM originaldm;
388: };
390: PETSC_EXTERN PetscErrorCode DMGetDMTS(DM,DMTS*);
391: PETSC_EXTERN PetscErrorCode DMGetDMTSWrite(DM,DMTS*);
392: PETSC_EXTERN PetscErrorCode DMCopyDMTS(DM,DM);
393: PETSC_EXTERN PetscErrorCode DMTSView(DMTS,PetscViewer);
394: PETSC_EXTERN PetscErrorCode DMTSLoad(DMTS,PetscViewer);
395: PETSC_EXTERN PetscErrorCode DMTSCopy(DMTS,DMTS);
397: typedef enum {TSEVENT_NONE,TSEVENT_LOCATED_INTERVAL,TSEVENT_PROCESSING,TSEVENT_ZERO,TSEVENT_RESET_NEXTSTEP} TSEventStatus;
399: struct _n_TSEvent {
400: PetscScalar *fvalue; /* value of event function at the end of the step*/
401: PetscScalar *fvalue_prev; /* value of event function at start of the step (left end-point of event interval) */
402: PetscReal ptime_prev; /* time at step start (left end-point of event interval) */
403: 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) */
404: PetscReal ptime_right; /* time on the right end-point of the event interval */
405: PetscScalar *fvalue_right; /* value of event function at the right end-point of the event interval */
406: PetscInt *side; /* Used for detecting repetition of end-point, -1 => left, +1 => right */
407: PetscReal timestep_prev; /* previous time step */
408: PetscReal timestep_posteventinterval; /* time step immediately after the event interval */
409: PetscBool *zerocrossing; /* Flag to signal zero crossing detection */
410: PetscErrorCode (*eventhandler)(TS,PetscReal,Vec,PetscScalar*,void*); /* User event handler function */
411: PetscErrorCode (*postevent)(TS,PetscInt,PetscInt[],PetscReal,Vec,PetscBool,void*); /* User post event function */
412: void *ctx; /* User context for event handler and post even functions */
413: PetscInt *direction; /* Zero crossing direction: 1 -> Going positive, -1 -> Going negative, 0 -> Any */
414: PetscBool *terminate; /* 1 -> Terminate time stepping, 0 -> continue */
415: PetscInt nevents; /* Number of events to handle */
416: PetscInt nevents_zero; /* Number of event zero detected */
417: PetscInt *events_zero; /* List of events that have reached zero */
418: PetscReal *vtol; /* Vector tolerances for event zero check */
419: TSEventStatus status; /* Event status */
420: PetscInt iterctr; /* Iteration counter */
421: PetscViewer monitor;
422: /* Struct to record the events */
423: struct {
424: PetscInt ctr; /* recorder counter */
425: PetscReal *time; /* Event times */
426: PetscInt *stepnum; /* Step numbers */
427: PetscInt *nevents; /* Number of events occuring at the event times */
428: PetscInt **eventidx; /* Local indices of the events in the event list */
429: } recorder;
430: PetscInt recsize; /* Size of recorder stack */
431: PetscInt refct; /* reference count */
432: };
434: PETSC_EXTERN PetscErrorCode TSEventInitialize(TSEvent,TS,PetscReal,Vec);
435: PETSC_EXTERN PetscErrorCode TSEventDestroy(TSEvent*);
436: PETSC_EXTERN PetscErrorCode TSEventHandler(TS);
437: PETSC_EXTERN PetscErrorCode TSAdjointEventHandler(TS);
439: PETSC_EXTERN PetscLogEvent TS_AdjointStep;
440: PETSC_EXTERN PetscLogEvent TS_Step;
441: PETSC_EXTERN PetscLogEvent TS_PseudoComputeTimeStep;
442: PETSC_EXTERN PetscLogEvent TS_FunctionEval;
443: PETSC_EXTERN PetscLogEvent TS_JacobianEval;
444: PETSC_EXTERN PetscLogEvent TS_ForwardStep;
446: typedef enum {TS_STEP_INCOMPLETE, /* vec_sol, ptime, etc point to beginning of step */
447: TS_STEP_PENDING, /* vec_sol advanced, but step has not been accepted yet */
448: TS_STEP_COMPLETE /* step accepted and ptime, steps, etc have been advanced */
449: } TSStepStatus;
451: struct _n_TSMonitorLGCtx {
452: PetscDrawLG lg;
453: PetscBool semilogy;
454: PetscInt howoften; /* when > 0 uses step % howoften, when negative only final solution plotted */
455: PetscInt ksp_its,snes_its;
456: char **names;
457: char **displaynames;
458: PetscInt ndisplayvariables;
459: PetscInt *displayvariables;
460: PetscReal *displayvalues;
461: PetscErrorCode (*transform)(void*,Vec,Vec*);
462: PetscErrorCode (*transformdestroy)(void*);
463: void *transformctx;
464: };
466: struct _n_TSMonitorSPCtx{
467: PetscDrawSP sp;
468: PetscInt howoften; /* when > 0 uses step % howoften, when negative only final solution plotted */
469: PetscInt ksp_its, snes_its;
470: };
472: struct _n_TSMonitorEnvelopeCtx {
473: Vec max,min;
474: };
476: /*
477: Checks if the user provide a TSSetIFunction() but an explicit method is called; generate an error in that case
478: */
479: PETSC_STATIC_INLINE PetscErrorCode TSCheckImplicitTerm(TS ts)
480: {
481: TSIFunction ifunction;
482: DM dm;
483: PetscErrorCode ierr;
486: TSGetDM(ts,&dm);
487: DMTSGetIFunction(dm,&ifunction,NULL);
488: if (ifunction) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_INCOMP,"You are attempting to use an explicit ODE integrator but provided an implicit function definition with TSSetIFunction()");
489: return(0);
490: }
492: PETSC_EXTERN PetscErrorCode TSGetRHSMats_Private(TS,Mat*,Mat*);
493: /* this is declared here as TSHistory is not public */
494: PETSC_EXTERN PetscErrorCode TSAdaptHistorySetTSHistory(TSAdapt,TSHistory,PetscBool);
496: PETSC_INTERN PetscErrorCode TSTrajectoryReconstruct_Private(TSTrajectory,TS,PetscReal,Vec,Vec);
497: PETSC_INTERN PetscErrorCode TSTrajectorySetUp_Basic(TSTrajectory,TS);
499: PETSC_EXTERN PetscLogEvent TSTrajectory_Set;
500: PETSC_EXTERN PetscLogEvent TSTrajectory_Get;
501: PETSC_EXTERN PetscLogEvent TSTrajectory_GetVecs;
502: PETSC_EXTERN PetscLogEvent TSTrajectory_DiskWrite;
503: PETSC_EXTERN PetscLogEvent TSTrajectory_DiskRead;
505: struct _n_TSMonitorDrawCtx {
506: PetscViewer viewer;
507: Vec initialsolution;
508: PetscBool showinitial;
509: PetscInt howoften; /* when > 0 uses step % howoften, when negative only final solution plotted */
510: PetscBool showtimestepandtime;
511: };
512: #endif