Actual source code: tsimpl.h

  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);
 30: PETSC_EXTERN PetscErrorCode TSIRKRegisterAll(void);

 32: typedef struct _TSOps *TSOps;

 34: struct _TSOps {
 35:   PetscErrorCode (*snesfunction)(SNES,Vec,Vec,TS);
 36:   PetscErrorCode (*snesjacobian)(SNES,Vec,Mat,Mat,TS);
 37:   PetscErrorCode (*setup)(TS);
 38:   PetscErrorCode (*step)(TS);
 39:   PetscErrorCode (*solve)(TS);
 40:   PetscErrorCode (*interpolate)(TS,PetscReal,Vec);
 41:   PetscErrorCode (*evaluatewlte)(TS,NormType,PetscInt*,PetscReal*);
 42:   PetscErrorCode (*evaluatestep)(TS,PetscInt,Vec,PetscBool*);
 43:   PetscErrorCode (*setfromoptions)(PetscOptionItems*,TS);
 44:   PetscErrorCode (*destroy)(TS);
 45:   PetscErrorCode (*view)(TS,PetscViewer);
 46:   PetscErrorCode (*reset)(TS);
 47:   PetscErrorCode (*linearstability)(TS,PetscReal,PetscReal,PetscReal*,PetscReal*);
 48:   PetscErrorCode (*load)(TS,PetscViewer);
 49:   PetscErrorCode (*rollback)(TS);
 50:   PetscErrorCode (*getstages)(TS,PetscInt*,Vec*[]);
 51:   PetscErrorCode (*adjointstep)(TS);
 52:   PetscErrorCode (*adjointsetup)(TS);
 53:   PetscErrorCode (*adjointreset)(TS);
 54:   PetscErrorCode (*adjointintegral)(TS);
 55:   PetscErrorCode (*forwardsetup)(TS);
 56:   PetscErrorCode (*forwardreset)(TS);
 57:   PetscErrorCode (*forwardstep)(TS);
 58:   PetscErrorCode (*forwardintegral)(TS);
 59:   PetscErrorCode (*forwardgetstages)(TS,PetscInt*,Mat*[]);
 60:   PetscErrorCode (*getsolutioncomponents)(TS,PetscInt*,Vec*);
 61:   PetscErrorCode (*getauxsolution)(TS,Vec*);
 62:   PetscErrorCode (*gettimeerror)(TS,PetscInt,Vec*);
 63:   PetscErrorCode (*settimeerror)(TS,Vec);
 64:   PetscErrorCode (*startingmethod) (TS);
 65:   PetscErrorCode (*initcondition)(TS,Vec);
 66:   PetscErrorCode (*exacterror)(TS,Vec,Vec);
 67: };

 69: /*
 70:    TSEvent - Abstract object to handle event monitoring
 71: */
 72: typedef struct _n_TSEvent *TSEvent;

 74: typedef struct _TSTrajectoryOps *TSTrajectoryOps;

 76: struct _TSTrajectoryOps {
 77:   PetscErrorCode (*view)(TSTrajectory,PetscViewer);
 78:   PetscErrorCode (*reset)(TSTrajectory);
 79:   PetscErrorCode (*destroy)(TSTrajectory);
 80:   PetscErrorCode (*set)(TSTrajectory,TS,PetscInt,PetscReal,Vec);
 81:   PetscErrorCode (*get)(TSTrajectory,TS,PetscInt,PetscReal*);
 82:   PetscErrorCode (*setfromoptions)(PetscOptionItems*,TSTrajectory);
 83:   PetscErrorCode (*setup)(TSTrajectory,TS);
 84: };

 86: /* TSHistory is an helper object that allows inquiring
 87:    the TSTrajectory by time and not by the step number only */
 88: typedef struct _n_TSHistory* TSHistory;

 90: struct _p_TSTrajectory {
 91:   PETSCHEADER(struct _TSTrajectoryOps);
 92:   TSHistory tsh;        /* associates times to unique step ids */
 93:   /* stores necessary data to reconstruct states and derivatives via Lagrangian interpolation */
 94:   struct {
 95:     PetscInt    order;  /* interpolation order */
 96:     Vec         *W;     /* work vectors */
 97:     PetscScalar *L;     /* workspace for Lagrange basis */
 98:     PetscReal   *T;     /* Lagrange times (stored) */
 99:     Vec         *WW;    /* just an array of pointers */
100:     PetscBool   *TT;    /* workspace for Lagrange */
101:     PetscReal   *TW;    /* Lagrange times (workspace) */

103:     /* caching */
104:     PetscBool caching;
105:     struct {
106:       PetscObjectId    id;
107:       PetscObjectState state;
108:       PetscReal        time;
109:       PetscInt         step;
110:     } Ucached;
111:     struct {
112:       PetscObjectId    id;
113:       PetscObjectState state;
114:       PetscReal        time;
115:       PetscInt         step;
116:     } Udotcached;
117:   } lag;
118:   Vec            U,Udot;                  /* used by TSTrajectory{Get|Restore}UpdatedHistoryVecs */
119:   PetscBool      usehistory;              /* whether to use TSHistory */
120:   PetscBool      solution_only;           /* whether we dump just the solution or also the stages */
121:   PetscBool      adjoint_solve_mode;      /* whether we will use the Trajectory inside a TSAdjointSolve() or not */
122:   PetscViewer    monitor;
123:   PetscInt       setupcalled;             /* true if setup has been called */
124:   PetscInt       recomps;                 /* counter for recomputations in the adjoint run */
125:   PetscInt       diskreads,diskwrites;    /* counters for disk checkpoint reads and writes */
126:   char           **names;                 /* the name of each variable; each process has only the local names */
127:   PetscBool      keepfiles;               /* keep the files generated during the run after the run is complete */
128:   char           *dirname,*filetemplate;  /* directory name and file name template for disk checkpoints */
129:   char           *dirfiletemplate;        /* complete directory and file name template for disk checkpoints */
130:   PetscErrorCode (*transform)(void*,Vec,Vec*);
131:   PetscErrorCode (*transformdestroy)(void*);
132:   void*          transformctx;
133:   void           *data;
134: };

136: typedef struct _TS_RHSSplitLink *TS_RHSSplitLink;
137: struct _TS_RHSSplitLink {
138:   TS              ts;
139:   char            *splitname;
140:   IS              is;
141:   TS_RHSSplitLink next;
142:   PetscLogEvent   event;
143: };

145: struct _p_TS {
146:   PETSCHEADER(struct _TSOps);
147:   TSProblemType  problem_type;
148:   TSEquationType equation_type;

150:   DM             dm;
151:   Vec            vec_sol; /* solution vector in first and second order equations */
152:   Vec            vec_dot; /* time derivative vector in second order equations */
153:   TSAdapt        adapt;
154:   TSAdaptType    default_adapt_type;
155:   TSEvent        event;

157:   /* ---------------- User (or PETSc) Provided stuff ---------------------*/
158:   PetscErrorCode (*monitor[MAXTSMONITORS])(TS,PetscInt,PetscReal,Vec,void*);
159:   PetscErrorCode (*monitordestroy[MAXTSMONITORS])(void**);
160:   void            *monitorcontext[MAXTSMONITORS];
161:   PetscInt         numbermonitors;
162:   PetscErrorCode (*adjointmonitor[MAXTSMONITORS])(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*);
163:   PetscErrorCode (*adjointmonitordestroy[MAXTSMONITORS])(void**);
164:   void            *adjointmonitorcontext[MAXTSMONITORS];
165:   PetscInt         numberadjointmonitors;
166:   PetscInt         monitorFrequency; /* Number of timesteps between monitor output */

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

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

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

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

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

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

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

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

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

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

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

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

276:   /* --------------------------------------------------------------------- */

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

288:   TSConvergedReason reason;
289:   PetscBool errorifstepfailed;
290:   PetscInt  reject,max_reject;
291:   TSExactFinalTimeOption exact_final_time;

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

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

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

308:   /* ---------------------- Quadrature integration support ---------------------------------*/
309:   TS quadraturets;
310: };

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

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

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

355:   TSIFunction ifunction;
356:   PetscErrorCode (*ifunctionview)(void*,PetscViewer);
357:   PetscErrorCode (*ifunctionload)(void**,PetscViewer);

359:   TSIJacobian ijacobian;
360:   PetscErrorCode (*ijacobianview)(void*,PetscViewer);
361:   PetscErrorCode (*ijacobianload)(void**,PetscViewer);

363:   TSI2Function i2function;
364:   TSI2Jacobian i2jacobian;

366:   TSTransientVariable transientvar;

368:   TSSolutionFunction solution;
369:   TSForcingFunction  forcing;

371:   PetscErrorCode (*destroy)(DMTS);
372:   PetscErrorCode (*duplicate)(DMTS,DMTS);
373: };

375: struct _p_DMTS {
376:   PETSCHEADER(struct _DMTSOps);
377:   void *rhsfunctionctx;
378:   void *rhsjacobianctx;

380:   void *ifunctionctx;
381:   void *ijacobianctx;

383:   void *i2functionctx;
384:   void *i2jacobianctx;

386:   void *transientvarctx;

388:   void *solutionctx;
389:   void *forcingctx;

391:   void *data;

393:   /* This is NOT reference counted. 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. Thus, if a user
395:    * only interacts directly with one level, e.g., using TSSetIFunction(), then coarse levels of a multilevel item
396:    * integrator are built, then the user changes the routine with another call to TSSetIFunction(), it automatically
397:    * propagates to all the levels. If instead, they get out a specific level and set the function on that level,
398:    * subsequent changes to the original level will no longer propagate to that level.
399:    */
400:   DM originaldm;
401: };

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

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

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

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

454: PETSC_EXTERN PetscLogEvent TS_AdjointStep;
455: PETSC_EXTERN PetscLogEvent TS_Step;
456: PETSC_EXTERN PetscLogEvent TS_PseudoComputeTimeStep;
457: PETSC_EXTERN PetscLogEvent TS_FunctionEval;
458: PETSC_EXTERN PetscLogEvent TS_JacobianEval;
459: PETSC_EXTERN PetscLogEvent TS_ForwardStep;

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

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

481: struct _n_TSMonitorSPCtx{
482:   PetscDrawSP sp;
483:   PetscInt    howoften; /* when > 0 uses step % howoften, when negative only final solution plotted */
484:   PetscInt    retain;   /* Retain n points plotted to show trajectories, or -1 for all points */
485:   PetscBool   phase;    /* Plot in phase space rather than coordinate space */
486:   PetscInt    ksp_its, snes_its;
487: };

489: struct _n_TSMonitorEnvelopeCtx {
490:   Vec max,min;
491: };

493: /*
494:     Checks if the user provide a TSSetIFunction() but an explicit method is called; generate an error in that case
495: */
496: static inline PetscErrorCode TSCheckImplicitTerm(TS ts)
497: {
498:   TSIFunction ifunction;
499:   DM          dm;

501:   TSGetDM(ts,&dm);
502:   DMTSGetIFunction(dm,&ifunction,NULL);
504:   return 0;
505: }

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

511: PETSC_INTERN PetscErrorCode TSTrajectoryReconstruct_Private(TSTrajectory,TS,PetscReal,Vec,Vec);
512: PETSC_INTERN PetscErrorCode TSTrajectorySetUp_Basic(TSTrajectory,TS);

514: PETSC_EXTERN PetscLogEvent TSTrajectory_Set;
515: PETSC_EXTERN PetscLogEvent TSTrajectory_Get;
516: PETSC_EXTERN PetscLogEvent TSTrajectory_GetVecs;
517: PETSC_EXTERN PetscLogEvent TSTrajectory_SetUp;
518: PETSC_EXTERN PetscLogEvent TSTrajectory_DiskWrite;
519: PETSC_EXTERN PetscLogEvent TSTrajectory_DiskRead;

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