Actual source code: tsimpl.h
petsc-master 2021-01-18
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: PetscErrorCode (*initcondition)(TS,Vec);
65: PetscErrorCode (*exacterror)(TS,Vec,Vec);
66: };
68: /*
69: TSEvent - Abstract object to handle event monitoring
70: */
71: typedef struct _n_TSEvent *TSEvent;
73: typedef struct _TSTrajectoryOps *TSTrajectoryOps;
75: struct _TSTrajectoryOps {
76: PetscErrorCode (*view)(TSTrajectory,PetscViewer);
77: PetscErrorCode (*reset)(TSTrajectory);
78: PetscErrorCode (*destroy)(TSTrajectory);
79: PetscErrorCode (*set)(TSTrajectory,TS,PetscInt,PetscReal,Vec);
80: PetscErrorCode (*get)(TSTrajectory,TS,PetscInt,PetscReal*);
81: PetscErrorCode (*setfromoptions)(PetscOptionItems*,TSTrajectory);
82: PetscErrorCode (*setup)(TSTrajectory,TS);
83: };
85: /* TSHistory is an helper object that allows inquiring
86: the TSTrajectory by time and not by the step number only */
87: typedef struct _n_TSHistory* TSHistory;
89: struct _p_TSTrajectory {
90: PETSCHEADER(struct _TSTrajectoryOps);
91: TSHistory tsh; /* associates times to unique step ids */
92: /* stores necessary data to reconstruct states and derivatives via Lagrangian interpolation */
93: struct {
94: PetscInt order; /* interpolation order */
95: Vec *W; /* work vectors */
96: PetscScalar *L; /* workspace for Lagrange basis */
97: PetscReal *T; /* Lagrange times (stored) */
98: Vec *WW; /* just an array of pointers */
99: PetscBool *TT; /* workspace for Lagrange */
100: PetscReal *TW; /* Lagrange times (workspace) */
102: /* caching */
103: PetscBool caching;
104: struct {
105: PetscObjectId id;
106: PetscObjectState state;
107: PetscReal time;
108: PetscInt step;
109: } Ucached;
110: struct {
111: PetscObjectId id;
112: PetscObjectState state;
113: PetscReal time;
114: PetscInt step;
115: } Udotcached;
116: } lag;
117: Vec U,Udot; /* used by TSTrajectory{Get|Restore}UpdatedHistoryVecs */
118: PetscBool usehistory; /* whether to use TSHistory */
119: PetscBool solution_only; /* whether we dump just the solution or also the stages */
120: PetscBool adjoint_solve_mode; /* whether we will use the Trajectory inside a TSAdjointSolve() or not */
121: PetscViewer monitor;
122: PetscInt setupcalled; /* true if setup has been called */
123: PetscInt recomps; /* counter for recomputations in the adjoint run */
124: PetscInt diskreads,diskwrites; /* counters for disk checkpoint reads and writes */
125: char **names; /* the name of each variable; each process has only the local names */
126: PetscBool keepfiles; /* keep the files generated during the run after the run is complete */
127: char *dirname,*filetemplate; /* directory name and file name template for disk checkpoints */
128: char *dirfiletemplate; /* complete directory and file name template for disk checkpoints */
129: PetscErrorCode (*transform)(void*,Vec,Vec*);
130: PetscErrorCode (*transformdestroy)(void*);
131: void* transformctx;
132: void *data;
133: };
135: typedef struct _TS_RHSSplitLink *TS_RHSSplitLink;
136: struct _TS_RHSSplitLink {
137: TS ts;
138: char *splitname;
139: IS is;
140: TS_RHSSplitLink next;
141: PetscLogEvent event;
142: };
144: struct _p_TS {
145: PETSCHEADER(struct _TSOps);
146: TSProblemType problem_type;
147: TSEquationType equation_type;
149: DM dm;
150: Vec vec_sol; /* solution vector in first and second order equations */
151: Vec vec_dot; /* time derivative vector in second order equations */
152: TSAdapt adapt;
153: TSAdaptType default_adapt_type;
154: TSEvent event;
156: /* ---------------- User (or PETSc) Provided stuff ---------------------*/
157: PetscErrorCode (*monitor[MAXTSMONITORS])(TS,PetscInt,PetscReal,Vec,void*);
158: PetscErrorCode (*monitordestroy[MAXTSMONITORS])(void**);
159: void *monitorcontext[MAXTSMONITORS];
160: PetscInt numbermonitors;
161: PetscErrorCode (*adjointmonitor[MAXTSMONITORS])(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*);
162: PetscErrorCode (*adjointmonitordestroy[MAXTSMONITORS])(void**);
163: void *adjointmonitorcontext[MAXTSMONITORS];
164: PetscInt numberadjointmonitors;
166: PetscErrorCode (*prestep)(TS);
167: PetscErrorCode (*prestage)(TS,PetscReal);
168: PetscErrorCode (*poststage)(TS,PetscReal,PetscInt,Vec*);
169: PetscErrorCode (*postevaluate)(TS);
170: PetscErrorCode (*poststep)(TS);
171: PetscErrorCode (*functiondomainerror)(TS,PetscReal,Vec,PetscBool*);
173: /* ---------------------- Sensitivity Analysis support -----------------*/
174: TSTrajectory trajectory; /* All solutions are kept here for the entire time integration process */
175: Vec *vecs_sensi; /* one vector for each cost function */
176: Vec *vecs_sensip;
177: PetscInt numcost; /* number of cost functions */
178: Vec vec_costintegral;
179: PetscInt adjointsetupcalled;
180: PetscInt adjoint_steps;
181: PetscInt adjoint_max_steps;
182: PetscBool adjoint_solve; /* immediately call TSAdjointSolve() after TSSolve() is complete */
183: PetscBool costintegralfwd; /* cost integral is evaluated in the forward run if true */
184: Vec vec_costintegrand; /* workspace for Adjoint computations */
185: Mat Jacp,Jacprhs;
186: void *ijacobianpctx,*rhsjacobianpctx;
187: void *costintegrandctx;
188: Vec *vecs_drdu;
189: Vec *vecs_drdp;
190: Vec vec_drdu_col,vec_drdp_col;
192: /* first-order adjoint */
193: PetscErrorCode (*rhsjacobianp)(TS,PetscReal,Vec,Mat,void*);
194: PetscErrorCode (*ijacobianp)(TS,PetscReal,Vec,Vec,PetscReal,Mat,void*);
195: PetscErrorCode (*costintegrand)(TS,PetscReal,Vec,Vec,void*);
196: PetscErrorCode (*drdufunction)(TS,PetscReal,Vec,Vec*,void*);
197: PetscErrorCode (*drdpfunction)(TS,PetscReal,Vec,Vec*,void*);
199: /* second-order adjoint */
200: Vec *vecs_sensi2;
201: Vec *vecs_sensi2p;
202: Vec vec_dir; /* directional vector for optimization */
203: Vec *vecs_fuu,*vecs_guu;
204: Vec *vecs_fup,*vecs_gup;
205: Vec *vecs_fpu,*vecs_gpu;
206: Vec *vecs_fpp,*vecs_gpp;
207: void *ihessianproductctx,*rhshessianproductctx;
208: PetscErrorCode (*ihessianproduct_fuu)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*);
209: PetscErrorCode (*ihessianproduct_fup)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*);
210: PetscErrorCode (*ihessianproduct_fpu)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*);
211: PetscErrorCode (*ihessianproduct_fpp)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*);
212: PetscErrorCode (*rhshessianproduct_guu)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*);
213: PetscErrorCode (*rhshessianproduct_gup)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*);
214: PetscErrorCode (*rhshessianproduct_gpu)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*);
215: PetscErrorCode (*rhshessianproduct_gpp)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*);
217: /* specific to forward sensitivity analysis */
218: Mat mat_sensip; /* matrix storing forward sensitivities */
219: Vec vec_sensip_col; /* space for a column of the sensip matrix */
220: Vec *vecs_integral_sensip; /* one vector for each integral */
221: PetscInt num_parameters;
222: PetscInt num_initialvalues;
223: void *vecsrhsjacobianpctx;
224: PetscInt forwardsetupcalled;
225: PetscBool forward_solve;
226: PetscErrorCode (*vecsrhsjacobianp)(TS,PetscReal,Vec,Vec*,void*);
228: /* ---------------------- IMEX support ---------------------------------*/
229: /* These extra slots are only used when the user provides both Implicit and RHS */
230: Mat Arhs; /* Right hand side matrix */
231: Mat Brhs; /* Right hand side preconditioning matrix */
232: Vec Frhs; /* Right hand side function value */
234: /* This is a general caching scheme to avoid recomputing the Jacobian at a place that has been previously been evaluated.
235: * The present use case is that TSComputeRHSFunctionLinear() evaluates the Jacobian once and we don't want it to be immeditely re-evaluated.
236: */
237: struct {
238: PetscReal time; /* The time at which the matrices were last evaluated */
239: PetscObjectId Xid; /* Unique ID of solution vector at which the Jacobian was last evaluated */
240: PetscObjectState Xstate; /* State of the solution vector */
241: MatStructure mstructure; /* The structure returned */
242: /* Flag to unshift Jacobian before calling the IJacobian or RHSJacobian functions. This is useful
243: * if the user would like to reuse (part of) the Jacobian from the last evaluation. */
244: PetscBool reuse;
245: PetscReal scale,shift;
246: } rhsjacobian;
248: struct {
249: PetscReal shift; /* The derivative of the lhs wrt to Xdot */
250: } ijacobian;
252: MatStructure axpy_pattern; /* information about the nonzero pattern of the RHS Jacobian in reference to the implicit Jacobian */
253: /* --------------------Nonlinear Iteration------------------------------*/
254: SNES snes;
255: PetscBool usessnes; /* Flag set by each TSType to indicate if the type actually uses a SNES;
256: this works around the design flaw that a SNES is ALWAYS created with TS even when it is not needed.*/
257: PetscInt ksp_its; /* total number of linear solver iterations */
258: PetscInt snes_its; /* total number of nonlinear solver iterations */
259: PetscInt num_snes_failures;
260: PetscInt max_snes_failures;
262: /* --- Logging --- */
263: PetscInt ifuncs,rhsfuncs,ijacs,rhsjacs;
265: /* --- Data that is unique to each particular solver --- */
266: PetscInt setupcalled; /* true if setup has been called */
267: void *data; /* implementationspecific data */
268: void *user; /* user context */
270: /* ------------------ Parameters -------------------------------------- */
271: PetscInt max_steps; /* max number of steps */
272: PetscReal max_time; /* max time allowed */
274: /* --------------------------------------------------------------------- */
276: PetscBool steprollback; /* flag to indicate that the step was rolled back */
277: PetscBool steprestart; /* flag to indicate that the timestepper has to discard any history and restart */
278: PetscInt steps; /* steps taken so far in all successive calls to TSSolve() */
279: PetscReal ptime; /* time at the start of the current step (stage time is internal if it exists) */
280: PetscReal time_step; /* current time increment */
281: PetscReal ptime_prev; /* time at the start of the previous step */
282: PetscReal ptime_prev_rollback; /* time at the start of the 2nd previous step to recover from rollback */
283: PetscReal solvetime; /* time at the conclusion of TSSolve() */
285: TSConvergedReason reason;
286: PetscBool errorifstepfailed;
287: PetscInt reject,max_reject;
288: TSExactFinalTimeOption exact_final_time;
290: PetscReal atol,rtol; /* Relative and absolute tolerance for local truncation error */
291: Vec vatol,vrtol; /* Relative and absolute tolerance in vector form */
292: PetscReal cfltime,cfltime_local;
294: PetscBool testjacobian;
295: PetscBool testjacobiantranspose;
296: /* ------------------- Default work-area management ------------------ */
297: PetscInt nwork;
298: Vec *work;
300: /* ---------------------- RHS splitting support ---------------------------------*/
301: PetscInt num_rhs_splits;
302: TS_RHSSplitLink tsrhssplit;
303: PetscBool use_splitrhsfunction;
305: /* ---------------------- Quadrature integration support ---------------------------------*/
306: TS quadraturets;
307: };
309: struct _TSAdaptOps {
310: PetscErrorCode (*choose)(TSAdapt,TS,PetscReal,PetscInt*,PetscReal*,PetscBool*,PetscReal*,PetscReal*,PetscReal*);
311: PetscErrorCode (*destroy)(TSAdapt);
312: PetscErrorCode (*reset)(TSAdapt);
313: PetscErrorCode (*view)(TSAdapt,PetscViewer);
314: PetscErrorCode (*setfromoptions)(PetscOptionItems*,TSAdapt);
315: PetscErrorCode (*load)(TSAdapt,PetscViewer);
316: };
318: struct _p_TSAdapt {
319: PETSCHEADER(struct _TSAdaptOps);
320: void *data;
321: PetscErrorCode (*checkstage)(TSAdapt,TS,PetscReal,Vec,PetscBool*);
322: struct {
323: PetscInt n; /* number of candidate schemes, including the one currently in use */
324: PetscBool inuse_set; /* the current scheme has been set */
325: const char *name[16]; /* name of the scheme */
326: PetscInt order[16]; /* classical order of each scheme */
327: PetscInt stageorder[16]; /* stage order of each scheme */
328: PetscReal ccfl[16]; /* stability limit relative to explicit Euler */
329: PetscReal cost[16]; /* relative measure of the amount of work required for each scheme */
330: } candidates;
331: PetscBool always_accept;
332: PetscReal safety; /* safety factor relative to target error/stability goal */
333: PetscReal reject_safety; /* extra safety factor if the last step was rejected */
334: PetscReal clip[2]; /* admissible time step decrease/increase factors */
335: PetscReal dt_min,dt_max; /* admissible minimum and maximum time step */
336: PetscReal ignore_max; /* minimum value of the solution to be considered by the adaptor */
337: PetscBool glee_use_local; /* GLEE adaptor uses global or local error */
338: PetscReal scale_solve_failed; /* scale step by this factor if solver (linear or nonlinear) fails. */
339: PetscReal matchstepfac[2]; /* factors to control the behaviour of matchstep */
340: NormType wnormtype;
341: PetscViewer monitor;
342: PetscInt timestepjustdecreased_delay; /* number of timesteps after a decrease in the timestep before the timestep can be increased */
343: PetscInt timestepjustdecreased;
344: };
346: typedef struct _p_DMTS *DMTS;
347: typedef struct _DMTSOps *DMTSOps;
348: struct _DMTSOps {
349: TSRHSFunction rhsfunction;
350: TSRHSJacobian rhsjacobian;
352: TSIFunction ifunction;
353: PetscErrorCode (*ifunctionview)(void*,PetscViewer);
354: PetscErrorCode (*ifunctionload)(void**,PetscViewer);
356: TSIJacobian ijacobian;
357: PetscErrorCode (*ijacobianview)(void*,PetscViewer);
358: PetscErrorCode (*ijacobianload)(void**,PetscViewer);
360: TSI2Function i2function;
361: TSI2Jacobian i2jacobian;
363: TSTransientVariable transientvar;
365: TSSolutionFunction solution;
366: TSForcingFunction forcing;
368: PetscErrorCode (*destroy)(DMTS);
369: PetscErrorCode (*duplicate)(DMTS,DMTS);
370: };
372: struct _p_DMTS {
373: PETSCHEADER(struct _DMTSOps);
374: void *rhsfunctionctx;
375: void *rhsjacobianctx;
377: void *ifunctionctx;
378: void *ijacobianctx;
380: void *i2functionctx;
381: void *i2jacobianctx;
383: void *transientvarctx;
385: void *solutionctx;
386: void *forcingctx;
388: void *data;
390: /* This is NOT reference counted. The DM on which this context was first created is cached here to implement one-way
391: * copy-on-write. When DMGetDMTSWrite() sees a request using a different DM, it makes a copy. Thus, if a user
392: * only interacts directly with one level, e.g., using TSSetIFunction(), then coarse levels of a multilevel item
393: * integrator are built, then the user changes the routine with another call to TSSetIFunction(), it automatically
394: * propagates to all the levels. If instead, they get out a specific level and set the function on that level,
395: * subsequent changes to the original level will no longer propagate to that level.
396: */
397: DM originaldm;
398: };
400: PETSC_EXTERN PetscErrorCode DMGetDMTS(DM,DMTS*);
401: PETSC_EXTERN PetscErrorCode DMGetDMTSWrite(DM,DMTS*);
402: PETSC_EXTERN PetscErrorCode DMCopyDMTS(DM,DM);
403: PETSC_EXTERN PetscErrorCode DMTSView(DMTS,PetscViewer);
404: PETSC_EXTERN PetscErrorCode DMTSLoad(DMTS,PetscViewer);
405: PETSC_EXTERN PetscErrorCode DMTSCopy(DMTS,DMTS);
407: typedef enum {TSEVENT_NONE,TSEVENT_LOCATED_INTERVAL,TSEVENT_PROCESSING,TSEVENT_ZERO,TSEVENT_RESET_NEXTSTEP} TSEventStatus;
409: struct _n_TSEvent {
410: PetscScalar *fvalue; /* value of event function at the end of the step*/
411: PetscScalar *fvalue_prev; /* value of event function at start of the step (left end-point of event interval) */
412: PetscReal ptime_prev; /* time at step start (left end-point of event interval) */
413: 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) */
414: PetscReal ptime_right; /* time on the right end-point of the event interval */
415: PetscScalar *fvalue_right; /* value of event function at the right end-point of the event interval */
416: PetscInt *side; /* Used for detecting repetition of end-point, -1 => left, +1 => right */
417: PetscReal timestep_prev; /* previous time step */
418: PetscReal timestep_posteventinterval; /* time step immediately after the event interval */
419: PetscBool *zerocrossing; /* Flag to signal zero crossing detection */
420: PetscErrorCode (*eventhandler)(TS,PetscReal,Vec,PetscScalar*,void*); /* User event handler function */
421: PetscErrorCode (*postevent)(TS,PetscInt,PetscInt[],PetscReal,Vec,PetscBool,void*); /* User post event function */
422: void *ctx; /* User context for event handler and post even functions */
423: PetscInt *direction; /* Zero crossing direction: 1 -> Going positive, -1 -> Going negative, 0 -> Any */
424: PetscBool *terminate; /* 1 -> Terminate time stepping, 0 -> continue */
425: PetscInt nevents; /* Number of events to handle */
426: PetscInt nevents_zero; /* Number of event zero detected */
427: PetscInt *events_zero; /* List of events that have reached zero */
428: PetscReal *vtol; /* Vector tolerances for event zero check */
429: TSEventStatus status; /* Event status */
430: PetscInt iterctr; /* Iteration counter */
431: PetscViewer monitor;
432: /* Struct to record the events */
433: struct {
434: PetscInt ctr; /* recorder counter */
435: PetscReal *time; /* Event times */
436: PetscInt *stepnum; /* Step numbers */
437: PetscInt *nevents; /* Number of events occuring at the event times */
438: PetscInt **eventidx; /* Local indices of the events in the event list */
439: } recorder;
440: PetscInt recsize; /* Size of recorder stack */
441: PetscInt refct; /* reference count */
442: };
444: PETSC_EXTERN PetscErrorCode TSEventInitialize(TSEvent,TS,PetscReal,Vec);
445: PETSC_EXTERN PetscErrorCode TSEventDestroy(TSEvent*);
446: PETSC_EXTERN PetscErrorCode TSEventHandler(TS);
447: PETSC_EXTERN PetscErrorCode TSAdjointEventHandler(TS);
449: PETSC_EXTERN PetscLogEvent TS_AdjointStep;
450: PETSC_EXTERN PetscLogEvent TS_Step;
451: PETSC_EXTERN PetscLogEvent TS_PseudoComputeTimeStep;
452: PETSC_EXTERN PetscLogEvent TS_FunctionEval;
453: PETSC_EXTERN PetscLogEvent TS_JacobianEval;
454: PETSC_EXTERN PetscLogEvent TS_ForwardStep;
456: typedef enum {TS_STEP_INCOMPLETE, /* vec_sol, ptime, etc point to beginning of step */
457: TS_STEP_PENDING, /* vec_sol advanced, but step has not been accepted yet */
458: TS_STEP_COMPLETE /* step accepted and ptime, steps, etc have been advanced */
459: } TSStepStatus;
461: struct _n_TSMonitorLGCtx {
462: PetscDrawLG lg;
463: PetscBool semilogy;
464: PetscInt howoften; /* when > 0 uses step % howoften, when negative only final solution plotted */
465: PetscInt ksp_its,snes_its;
466: char **names;
467: char **displaynames;
468: PetscInt ndisplayvariables;
469: PetscInt *displayvariables;
470: PetscReal *displayvalues;
471: PetscErrorCode (*transform)(void*,Vec,Vec*);
472: PetscErrorCode (*transformdestroy)(void*);
473: void *transformctx;
474: };
476: struct _n_TSMonitorSPCtx{
477: PetscDrawSP sp;
478: PetscInt howoften; /* when > 0 uses step % howoften, when negative only final solution plotted */
479: PetscInt ksp_its, snes_its;
480: };
482: struct _n_TSMonitorEnvelopeCtx {
483: Vec max,min;
484: };
486: /*
487: Checks if the user provide a TSSetIFunction() but an explicit method is called; generate an error in that case
488: */
489: PETSC_STATIC_INLINE PetscErrorCode TSCheckImplicitTerm(TS ts)
490: {
491: TSIFunction ifunction;
492: DM dm;
493: PetscErrorCode ierr;
496: TSGetDM(ts,&dm);
497: DMTSGetIFunction(dm,&ifunction,NULL);
498: 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()");
499: return(0);
500: }
502: PETSC_EXTERN PetscErrorCode TSGetRHSMats_Private(TS,Mat*,Mat*);
503: /* this is declared here as TSHistory is not public */
504: PETSC_EXTERN PetscErrorCode TSAdaptHistorySetTSHistory(TSAdapt,TSHistory,PetscBool);
506: PETSC_INTERN PetscErrorCode TSTrajectoryReconstruct_Private(TSTrajectory,TS,PetscReal,Vec,Vec);
507: PETSC_INTERN PetscErrorCode TSTrajectorySetUp_Basic(TSTrajectory,TS);
509: PETSC_EXTERN PetscLogEvent TSTrajectory_Set;
510: PETSC_EXTERN PetscLogEvent TSTrajectory_Get;
511: PETSC_EXTERN PetscLogEvent TSTrajectory_GetVecs;
512: PETSC_EXTERN PetscLogEvent TSTrajectory_DiskWrite;
513: PETSC_EXTERN PetscLogEvent TSTrajectory_DiskRead;
515: struct _n_TSMonitorDrawCtx {
516: PetscViewer viewer;
517: Vec initialsolution;
518: PetscBool showinitial;
519: PetscInt howoften; /* when > 0 uses step % howoften, when negative only final solution plotted */
520: PetscBool showtimestepandtime;
521: };
522: #endif