Actual source code: tsimpl.h

petsc-main 2021-04-20
Report Typos and Errors
  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;
165:   PetscInt         monitorFrequency; /* Number of timesteps between monitor output */

167:   PetscErrorCode (*prestep)(TS);
168:   PetscErrorCode (*prestage)(TS,PetscReal);
169:   PetscErrorCode (*poststage)(TS,PetscReal,PetscInt,Vec*);
170:   PetscErrorCode (*postevaluate)(TS);
171:   PetscErrorCode (*poststep)(TS);
172:   PetscErrorCode (*functiondomainerror)(TS,PetscReal,Vec,PetscBool*);

174:   /* ---------------------- Sensitivity Analysis support -----------------*/
175:   TSTrajectory trajectory;          /* All solutions are kept here for the entire time integration process */
176:   Vec       *vecs_sensi;            /* one vector for each cost function */
177:   Vec       *vecs_sensip;
178:   PetscInt  numcost;                /* number of cost functions */
179:   Vec       vec_costintegral;
180:   PetscInt  adjointsetupcalled;
181:   PetscInt  adjoint_steps;
182:   PetscInt  adjoint_max_steps;
183:   PetscBool adjoint_solve;          /* immediately call TSAdjointSolve() after TSSolve() is complete */
184:   PetscBool costintegralfwd;        /* cost integral is evaluated in the forward run if true */
185:   Vec       vec_costintegrand;      /* workspace for Adjoint computations */
186:   Mat       Jacp,Jacprhs;
187:   void      *ijacobianpctx,*rhsjacobianpctx;
188:   void      *costintegrandctx;
189:   Vec       *vecs_drdu;
190:   Vec       *vecs_drdp;
191:   Vec       vec_drdu_col,vec_drdp_col;

193:   /* first-order adjoint */
194:   PetscErrorCode (*rhsjacobianp)(TS,PetscReal,Vec,Mat,void*);
195:   PetscErrorCode (*ijacobianp)(TS,PetscReal,Vec,Vec,PetscReal,Mat,void*);
196:   PetscErrorCode (*costintegrand)(TS,PetscReal,Vec,Vec,void*);
197:   PetscErrorCode (*drdufunction)(TS,PetscReal,Vec,Vec*,void*);
198:   PetscErrorCode (*drdpfunction)(TS,PetscReal,Vec,Vec*,void*);

200:   /* second-order adjoint */
201:   Vec *vecs_sensi2;
202:   Vec *vecs_sensi2p;
203:   Vec vec_dir; /* directional vector for optimization */
204:   Vec *vecs_fuu,*vecs_guu;
205:   Vec *vecs_fup,*vecs_gup;
206:   Vec *vecs_fpu,*vecs_gpu;
207:   Vec *vecs_fpp,*vecs_gpp;
208:   void *ihessianproductctx,*rhshessianproductctx;
209:   PetscErrorCode (*ihessianproduct_fuu)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*);
210:   PetscErrorCode (*ihessianproduct_fup)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*);
211:   PetscErrorCode (*ihessianproduct_fpu)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*);
212:   PetscErrorCode (*ihessianproduct_fpp)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*);
213:   PetscErrorCode (*rhshessianproduct_guu)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*);
214:   PetscErrorCode (*rhshessianproduct_gup)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*);
215:   PetscErrorCode (*rhshessianproduct_gpu)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*);
216:   PetscErrorCode (*rhshessianproduct_gpp)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*);

218:   /* specific to forward sensitivity analysis */
219:   Mat       mat_sensip;              /* matrix storing forward sensitivities */
220:   Vec       vec_sensip_col;          /* space for a column of the sensip matrix */
221:   Vec       *vecs_integral_sensip;   /* one vector for each integral */
222:   PetscInt  num_parameters;
223:   PetscInt  num_initialvalues;
224:   void      *vecsrhsjacobianpctx;
225:   PetscInt  forwardsetupcalled;
226:   PetscBool forward_solve;
227:   PetscErrorCode (*vecsrhsjacobianp)(TS,PetscReal,Vec,Vec*,void*);

229:   /* ---------------------- IMEX support ---------------------------------*/
230:   /* These extra slots are only used when the user provides both Implicit and RHS */
231:   Mat Arhs;     /* Right hand side matrix */
232:   Mat Brhs;     /* Right hand side preconditioning matrix */
233:   Vec Frhs;     /* Right hand side function value */

235:   /* This is a general caching scheme to avoid recomputing the Jacobian at a place that has been previously been evaluated.
236:    * The present use case is that TSComputeRHSFunctionLinear() evaluates the Jacobian once and we don't want it to be immeditely re-evaluated.
237:    */
238:   struct {
239:     PetscReal        time;          /* The time at which the matrices were last evaluated */
240:     PetscObjectId    Xid;           /* Unique ID of solution vector at which the Jacobian was last evaluated */
241:     PetscObjectState Xstate;        /* State of the solution vector */
242:     MatStructure     mstructure;    /* The structure returned */
243:     /* Flag to unshift Jacobian before calling the IJacobian or RHSJacobian functions.  This is useful
244:      * if the user would like to reuse (part of) the Jacobian from the last evaluation. */
245:     PetscBool        reuse;
246:     PetscReal        scale,shift;
247:   } rhsjacobian;

249:   struct {
250:     PetscReal shift;            /* The derivative of the lhs wrt to Xdot */
251:   } ijacobian;

253:   MatStructure  axpy_pattern;  /* information about the nonzero pattern of the RHS Jacobian in reference to the implicit Jacobian */
254:   /* --------------------Nonlinear Iteration------------------------------*/
255:   SNES     snes;
256:   PetscBool usessnes;   /* Flag set by each TSType to indicate if the type actually uses a SNES;
257:                            this works around the design flaw that a SNES is ALWAYS created with TS even when it is not needed.*/
258:   PetscInt ksp_its;                /* total number of linear solver iterations */
259:   PetscInt snes_its;               /* total number of nonlinear solver iterations */
260:   PetscInt num_snes_failures;
261:   PetscInt max_snes_failures;

263:   /* --- Logging --- */
264:   PetscInt ifuncs,rhsfuncs,ijacs,rhsjacs;

266:   /* --- Data that is unique to each particular solver --- */
267:   PetscInt setupcalled;             /* true if setup has been called */
268:   void     *data;                   /* implementationspecific data */
269:   void     *user;                   /* user context */

271:   /* ------------------  Parameters -------------------------------------- */
272:   PetscInt  max_steps;              /* max number of steps */
273:   PetscReal max_time;               /* max time allowed */

275:   /* --------------------------------------------------------------------- */

277:   PetscBool steprollback;           /* flag to indicate that the step was rolled back */
278:   PetscBool steprestart;            /* flag to indicate that the timestepper has to discard any history and restart */
279:   PetscInt  steps;                  /* steps taken so far in all successive calls to TSSolve() */
280:   PetscReal ptime;                  /* time at the start of the current step (stage time is internal if it exists) */
281:   PetscReal time_step;              /* current time increment */
282:   PetscReal ptime_prev;             /* time at the start of the previous step */
283:   PetscReal ptime_prev_rollback;    /* time at the start of the 2nd previous step to recover from rollback */
284:   PetscReal solvetime;              /* time at the conclusion of TSSolve() */

286:   TSConvergedReason reason;
287:   PetscBool errorifstepfailed;
288:   PetscInt  reject,max_reject;
289:   TSExactFinalTimeOption exact_final_time;

291:   PetscReal atol,rtol;              /* Relative and absolute tolerance for local truncation error */
292:   Vec       vatol,vrtol;            /* Relative and absolute tolerance in vector form */
293:   PetscReal cfltime,cfltime_local;

295:   PetscBool testjacobian;
296:   PetscBool testjacobiantranspose;
297:   /* ------------------- Default work-area management ------------------ */
298:   PetscInt nwork;
299:   Vec      *work;

301:   /* ---------------------- RHS splitting support ---------------------------------*/
302:   PetscInt        num_rhs_splits;
303:   TS_RHSSplitLink tsrhssplit;
304:   PetscBool       use_splitrhsfunction;

306:   /* ---------------------- Quadrature integration support ---------------------------------*/
307:   TS quadraturets;
308: };

310: struct _TSAdaptOps {
311:   PetscErrorCode (*choose)(TSAdapt,TS,PetscReal,PetscInt*,PetscReal*,PetscBool*,PetscReal*,PetscReal*,PetscReal*);
312:   PetscErrorCode (*destroy)(TSAdapt);
313:   PetscErrorCode (*reset)(TSAdapt);
314:   PetscErrorCode (*view)(TSAdapt,PetscViewer);
315:   PetscErrorCode (*setfromoptions)(PetscOptionItems*,TSAdapt);
316:   PetscErrorCode (*load)(TSAdapt,PetscViewer);
317: };

319: struct _p_TSAdapt {
320:   PETSCHEADER(struct _TSAdaptOps);
321:   void *data;
322:   PetscErrorCode (*checkstage)(TSAdapt,TS,PetscReal,Vec,PetscBool*);
323:   struct {
324:     PetscInt   n;                /* number of candidate schemes, including the one currently in use */
325:     PetscBool  inuse_set;        /* the current scheme has been set */
326:     const char *name[16];        /* name of the scheme */
327:     PetscInt   order[16];        /* classical order of each scheme */
328:     PetscInt   stageorder[16];   /* stage order of each scheme */
329:     PetscReal  ccfl[16];         /* stability limit relative to explicit Euler */
330:     PetscReal  cost[16];         /* relative measure of the amount of work required for each scheme */
331:   } candidates;
332:   PetscBool   always_accept;
333:   PetscReal   safety;             /* safety factor relative to target error/stability goal */
334:   PetscReal   reject_safety;      /* extra safety factor if the last step was rejected */
335:   PetscReal   clip[2];            /* admissible time step decrease/increase factors */
336:   PetscReal   dt_min,dt_max;      /* admissible minimum and maximum time step */
337:   PetscReal   ignore_max;         /* minimum value of the solution to be considered by the adaptor */
338:   PetscBool   glee_use_local;     /* GLEE adaptor uses global or local error */
339:   PetscReal   scale_solve_failed; /* scale step by this factor if solver (linear or nonlinear) fails. */
340:   PetscReal   matchstepfac[2];    /* factors to control the behaviour of matchstep */
341:   NormType    wnormtype;
342:   PetscViewer monitor;
343:   PetscInt    timestepjustdecreased_delay; /* number of timesteps after a decrease in the timestep before the timestep can be increased */
344:   PetscInt    timestepjustdecreased;
345: };

347: typedef struct _p_DMTS *DMTS;
348: typedef struct _DMTSOps *DMTSOps;
349: struct _DMTSOps {
350:   TSRHSFunction rhsfunction;
351:   TSRHSJacobian rhsjacobian;

353:   TSIFunction ifunction;
354:   PetscErrorCode (*ifunctionview)(void*,PetscViewer);
355:   PetscErrorCode (*ifunctionload)(void**,PetscViewer);

357:   TSIJacobian ijacobian;
358:   PetscErrorCode (*ijacobianview)(void*,PetscViewer);
359:   PetscErrorCode (*ijacobianload)(void**,PetscViewer);

361:   TSI2Function i2function;
362:   TSI2Jacobian i2jacobian;

364:   TSTransientVariable transientvar;

366:   TSSolutionFunction solution;
367:   TSForcingFunction  forcing;

369:   PetscErrorCode (*destroy)(DMTS);
370:   PetscErrorCode (*duplicate)(DMTS,DMTS);
371: };

373: struct _p_DMTS {
374:   PETSCHEADER(struct _DMTSOps);
375:   void *rhsfunctionctx;
376:   void *rhsjacobianctx;

378:   void *ifunctionctx;
379:   void *ijacobianctx;

381:   void *i2functionctx;
382:   void *i2jacobianctx;

384:   void *transientvarctx;

386:   void *solutionctx;
387:   void *forcingctx;

389:   void *data;

391:   /* This is NOT reference counted. The DM on which this context was first created is cached here to implement one-way
392:    * copy-on-write. When DMGetDMTSWrite() sees a request using a different DM, it makes a copy. Thus, if a user
393:    * only interacts directly with one level, e.g., using TSSetIFunction(), then coarse levels of a multilevel item
394:    * integrator are built, then the user changes the routine with another call to TSSetIFunction(), it automatically
395:    * propagates to all the levels. If instead, they get out a specific level and set the function on that level,
396:    * subsequent changes to the original level will no longer propagate to that level.
397:    */
398:   DM originaldm;
399: };

401: PETSC_EXTERN PetscErrorCode DMGetDMTS(DM,DMTS*);
402: PETSC_EXTERN PetscErrorCode DMGetDMTSWrite(DM,DMTS*);
403: PETSC_EXTERN PetscErrorCode DMCopyDMTS(DM,DM);
404: PETSC_EXTERN PetscErrorCode DMTSView(DMTS,PetscViewer);
405: PETSC_EXTERN PetscErrorCode DMTSLoad(DMTS,PetscViewer);
406: PETSC_EXTERN PetscErrorCode DMTSCopy(DMTS,DMTS);

408: typedef enum {TSEVENT_NONE,TSEVENT_LOCATED_INTERVAL,TSEVENT_PROCESSING,TSEVENT_ZERO,TSEVENT_RESET_NEXTSTEP} TSEventStatus;

410: struct _n_TSEvent {
411:   PetscScalar    *fvalue;          /* value of event function at the end of the step*/
412:   PetscScalar    *fvalue_prev;     /* value of event function at start of the step (left end-point of event interval) */
413:   PetscReal       ptime_prev;      /* time at step start (left end-point of event interval) */
414:   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) */
415:   PetscReal       ptime_right;     /* time on the right end-point of the event interval */
416:   PetscScalar    *fvalue_right;    /* value of event function at the right end-point of the event interval */
417:   PetscInt       *side;            /* Used for detecting repetition of end-point, -1 => left, +1 => right */
418:   PetscReal       timestep_prev;   /* previous time step */
419:   PetscReal       timestep_posteventinterval;  /* time step immediately after the event interval */
420:   PetscReal       timestep_postevent;  /* time step immediately after the event */
421:   PetscReal       timestep_min;    /* Minimum time step */
422:   PetscBool      *zerocrossing;    /* Flag to signal zero crossing detection */
423:   PetscErrorCode  (*eventhandler)(TS,PetscReal,Vec,PetscScalar*,void*); /* User event handler function */
424:   PetscErrorCode  (*postevent)(TS,PetscInt,PetscInt[],PetscReal,Vec,PetscBool,void*); /* User post event function */
425:   void           *ctx;              /* User context for event handler and post even functions */
426:   PetscInt       *direction;        /* Zero crossing direction: 1 -> Going positive, -1 -> Going negative, 0 -> Any */
427:   PetscBool      *terminate;        /* 1 -> Terminate time stepping, 0 -> continue */
428:   PetscInt        nevents;          /* Number of events to handle */
429:   PetscInt        nevents_zero;     /* Number of event zero detected */
430:   PetscInt       *events_zero;      /* List of events that have reached zero */
431:   PetscReal      *vtol;             /* Vector tolerances for event zero check */
432:   TSEventStatus   status;           /* Event status */
433:   PetscInt        iterctr;          /* Iteration counter */
434:   PetscViewer     monitor;
435:   /* Struct to record the events */
436:   struct {
437:     PetscInt  ctr;        /* recorder counter */
438:     PetscReal *time;      /* Event times */
439:     PetscInt  *stepnum;   /* Step numbers */
440:     PetscInt  *nevents;   /* Number of events occuring at the event times */
441:     PetscInt  **eventidx; /* Local indices of the events in the event list */
442:   } recorder;
443:   PetscInt  recsize; /* Size of recorder stack */
444:   PetscInt  refct; /* reference count */
445: };

447: PETSC_EXTERN PetscErrorCode TSEventInitialize(TSEvent,TS,PetscReal,Vec);
448: PETSC_EXTERN PetscErrorCode TSEventDestroy(TSEvent*);
449: PETSC_EXTERN PetscErrorCode TSEventHandler(TS);
450: PETSC_EXTERN PetscErrorCode TSAdjointEventHandler(TS);

452: PETSC_EXTERN PetscLogEvent TS_AdjointStep;
453: PETSC_EXTERN PetscLogEvent TS_Step;
454: PETSC_EXTERN PetscLogEvent TS_PseudoComputeTimeStep;
455: PETSC_EXTERN PetscLogEvent TS_FunctionEval;
456: PETSC_EXTERN PetscLogEvent TS_JacobianEval;
457: PETSC_EXTERN PetscLogEvent TS_ForwardStep;

459: typedef enum {TS_STEP_INCOMPLETE, /* vec_sol, ptime, etc point to beginning of step */
460:               TS_STEP_PENDING,    /* vec_sol advanced, but step has not been accepted yet */
461:               TS_STEP_COMPLETE    /* step accepted and ptime, steps, etc have been advanced */
462: } TSStepStatus;

464: struct _n_TSMonitorLGCtx {
465:   PetscDrawLG    lg;
466:   PetscBool      semilogy;
467:   PetscInt       howoften;  /* when > 0 uses step % howoften, when negative only final solution plotted */
468:   PetscInt       ksp_its,snes_its;
469:   char           **names;
470:   char           **displaynames;
471:   PetscInt       ndisplayvariables;
472:   PetscInt       *displayvariables;
473:   PetscReal      *displayvalues;
474:   PetscErrorCode (*transform)(void*,Vec,Vec*);
475:   PetscErrorCode (*transformdestroy)(void*);
476:   void           *transformctx;
477: };

479: struct _n_TSMonitorSPCtx{
480:   PetscDrawSP    sp;
481:   PetscInt       howoften; /* when > 0 uses step % howoften, when negative only final solution plotted */
482:   PetscInt       ksp_its, snes_its;
483: };

485: struct _n_TSMonitorEnvelopeCtx {
486:   Vec max,min;
487: };

489: /*
490:     Checks if the user provide a TSSetIFunction() but an explicit method is called; generate an error in that case
491: */
492: PETSC_STATIC_INLINE PetscErrorCode TSCheckImplicitTerm(TS ts)
493: {
494:   TSIFunction      ifunction;
495:   DM               dm;
496:   PetscErrorCode   ierr;

499:   TSGetDM(ts,&dm);
500:   DMTSGetIFunction(dm,&ifunction,NULL);
501:   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()");
502:   return(0);
503: }

505: PETSC_EXTERN PetscErrorCode TSGetRHSMats_Private(TS,Mat*,Mat*);
506: /* this is declared here as TSHistory is not public */
507: PETSC_EXTERN PetscErrorCode TSAdaptHistorySetTSHistory(TSAdapt,TSHistory,PetscBool);

509: PETSC_INTERN PetscErrorCode TSTrajectoryReconstruct_Private(TSTrajectory,TS,PetscReal,Vec,Vec);
510: PETSC_INTERN PetscErrorCode TSTrajectorySetUp_Basic(TSTrajectory,TS);

512: PETSC_EXTERN PetscLogEvent TSTrajectory_Set;
513: PETSC_EXTERN PetscLogEvent TSTrajectory_Get;
514: PETSC_EXTERN PetscLogEvent TSTrajectory_GetVecs;
515: PETSC_EXTERN PetscLogEvent TSTrajectory_DiskWrite;
516: PETSC_EXTERN PetscLogEvent TSTrajectory_DiskRead;

518: struct _n_TSMonitorDrawCtx {
519:   PetscViewer   viewer;
520:   Vec           initialsolution;
521:   PetscBool     showinitial;
522:   PetscInt      howoften;  /* when > 0 uses step % howoften, when negative only final solution plotted */
523:   PetscBool     showtimestepandtime;
524: };
525: #endif