Actual source code: tsimpl.h
petsc-3.7.1 2016-05-15
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 TSARKIMEXRegisterAll(void);
26: PETSC_EXTERN PetscErrorCode TSRosWRegisterAll(void);
27: PETSC_EXTERN PetscErrorCode TSGLRegisterAll(void);
28: PETSC_EXTERN PetscErrorCode TSGLAdaptRegisterAll(void);
30: typedef struct _TSOps *TSOps;
32: struct _TSOps {
33: PetscErrorCode (*snesfunction)(SNES,Vec,Vec,TS);
34: PetscErrorCode (*snesjacobian)(SNES,Vec,Mat,Mat,TS);
35: PetscErrorCode (*setup)(TS);
36: PetscErrorCode (*step)(TS);
37: PetscErrorCode (*solve)(TS);
38: PetscErrorCode (*interpolate)(TS,PetscReal,Vec);
39: PetscErrorCode (*evaluatewlte)(TS,NormType,PetscInt*,PetscReal*);
40: PetscErrorCode (*evaluatestep)(TS,PetscInt,Vec,PetscBool*);
41: PetscErrorCode (*setfromoptions)(PetscOptionItems*,TS);
42: PetscErrorCode (*destroy)(TS);
43: PetscErrorCode (*view)(TS,PetscViewer);
44: PetscErrorCode (*reset)(TS);
45: PetscErrorCode (*linearstability)(TS,PetscReal,PetscReal,PetscReal*,PetscReal*);
46: PetscErrorCode (*load)(TS,PetscViewer);
47: PetscErrorCode (*rollback)(TS);
48: PetscErrorCode (*getstages)(TS,PetscInt*,Vec**);
49: PetscErrorCode (*adjointstep)(TS);
50: PetscErrorCode (*adjointsetup)(TS);
51: PetscErrorCode (*adjointintegral)(TS);
52: PetscErrorCode (*forwardintegral)(TS);
53: };
55: /*
56: TSEvent - Abstract object to handle event monitoring
57: */
58: typedef struct _n_TSEvent *TSEvent;
60: typedef struct _TSTrajectoryOps *TSTrajectoryOps;
62: struct _TSTrajectoryOps {
63: PetscErrorCode (*view)(TSTrajectory,PetscViewer);
64: PetscErrorCode (*destroy)(TSTrajectory);
65: PetscErrorCode (*set)(TSTrajectory,TS,PetscInt,PetscReal,Vec);
66: PetscErrorCode (*get)(TSTrajectory,TS,PetscInt,PetscReal*);
67: PetscErrorCode (*setfromoptions)(PetscOptionItems*,TSTrajectory);
68: PetscErrorCode (*setup)(TSTrajectory,TS);
69: };
71: struct _p_TSTrajectory {
72: PETSCHEADER(struct _TSTrajectoryOps);
73: PetscInt setupcalled; /* true if setup has been called */
74: PetscInt recomps; /* counter for recomputations in the adjoint run */
75: PetscInt diskreads,diskwrites; /* counters for disk checkpoint reads and writes */
76: void *data;
77: };
79: struct _p_TS {
80: PETSCHEADER(struct _TSOps);
81: TSProblemType problem_type;
82: TSEquationType equation_type;
84: DM dm;
85: Vec vec_sol; /* solution vector in first and second order equations */
86: Vec vec_dot; /* time derivative vector in second order equations */
87: TSAdapt adapt;
88: TSEvent event;
90: /* ---------------- User (or PETSc) Provided stuff ---------------------*/
91: PetscErrorCode (*monitor[MAXTSMONITORS])(TS,PetscInt,PetscReal,Vec,void*);
92: PetscErrorCode (*monitordestroy[MAXTSMONITORS])(void**);
93: void *monitorcontext[MAXTSMONITORS];
94: PetscInt numbermonitors;
95: PetscErrorCode (*adjointmonitor[MAXTSMONITORS])(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*);
96: PetscErrorCode (*adjointmonitordestroy[MAXTSMONITORS])(void**);
97: void *adjointmonitorcontext[MAXTSMONITORS];
98: PetscInt numberadjointmonitors;
100: PetscErrorCode (*prestep)(TS);
101: PetscErrorCode (*prestage)(TS,PetscReal);
102: PetscErrorCode (*poststage)(TS,PetscReal,PetscInt,Vec*);
103: PetscErrorCode (*poststep)(TS);
104: PetscErrorCode (*functiondomainerror)(TS,PetscReal,Vec,PetscBool*);
106: /* ---------------------- Sensitivity Analysis support -----------------*/
107: TSTrajectory trajectory; /* All solutions are kept here for the entire time integration process */
108: Vec *vecs_sensi; /* one vector for each cost function */
109: Vec *vecs_sensip;
110: PetscInt numcost; /* number of cost functions */
111: Vec vec_costintegral;
112: PetscInt adjointsetupcalled;
113: PetscInt adjoint_max_steps;
114: PetscBool adjoint_solve; /* immediately call TSAdjointSolve() after TSSolve() is complete */
115: PetscBool costintegralfwd; /* cost integral is evaluated in the forward run if true */
116: Vec vec_costintegrand; /* workspace for Adjoint computations */
117: Mat Jacp;
118: void *rhsjacobianpctx;
119: void *costintegrandctx;
120: Vec *vecs_drdy;
121: Vec *vecs_drdp;
123: PetscErrorCode (*rhsjacobianp)(TS,PetscReal,Vec,Mat,void*);
124: PetscErrorCode (*costintegrand)(TS,PetscReal,Vec,Vec,void*);
125: PetscErrorCode (*drdyfunction)(TS,PetscReal,Vec,Vec*,void*);
126: PetscErrorCode (*drdpfunction)(TS,PetscReal,Vec,Vec*,void*);
128: /* ---------------------- IMEX support ---------------------------------*/
129: /* These extra slots are only used when the user provides both Implicit and RHS */
130: Mat Arhs; /* Right hand side matrix */
131: Mat Brhs; /* Right hand side preconditioning matrix */
132: Vec Frhs; /* Right hand side function value */
134: /* This is a general caching scheme to avoid recomputing the Jacobian at a place that has been previously been evaluated.
135: * The present use case is that TSComputeRHSFunctionLinear() evaluates the Jacobian once and we don't want it to be immeditely re-evaluated.
136: */
137: struct {
138: PetscReal time; /* The time at which the matrices were last evaluated */
139: Vec X; /* Solution vector at which the Jacobian was last evaluated */
140: PetscObjectState Xstate; /* State of the solution vector */
141: MatStructure mstructure; /* The structure returned */
142: /* Flag to unshift Jacobian before calling the IJacobian or RHSJacobian functions. This is useful
143: * if the user would like to reuse (part of) the Jacobian from the last evaluation. */
144: PetscBool reuse;
145: PetscReal scale,shift;
146: } rhsjacobian;
148: struct {
149: PetscReal shift; /* The derivative of the lhs wrt to Xdot */
150: } ijacobian;
152: /* --------------------Nonlinear Iteration------------------------------*/
153: SNES snes;
154: PetscInt ksp_its; /* total number of linear solver iterations */
155: PetscInt snes_its; /* total number of nonlinear solver iterations */
156: PetscInt num_snes_failures;
157: PetscInt max_snes_failures;
159: /* --- Data that is unique to each particular solver --- */
160: PetscInt setupcalled; /* true if setup has been called */
161: void *data; /* implementationspecific data */
162: void *user; /* user context */
164: /* ------------------ Parameters -------------------------------------- */
165: PetscInt max_steps; /* max number of steps */
166: PetscReal max_time; /* max time allowed */
168: /* --------------------------------------------------------------------- */
170: PetscBool steprollback; /* flag to indicate that the step was rolled back */
171: PetscBool steprestart; /* flag to indicate that the timestepper has to discard any history and restart */
172: PetscInt steps; /* steps taken so far in latest call to TSSolve() */
173: PetscInt total_steps; /* steps taken in all calls to TSSolve() since the TS was created or since TSSetUp() was called */
174: PetscReal ptime; /* time at the start of the current step (stage time is internal if it exists) */
175: PetscReal time_step; /* current time increment */
176: PetscReal ptime_prev; /* time at the start of the previous step */
177: PetscReal ptime_prev_rollback; /* time at the start of the 2nd previous step to recover from rollback */
178: PetscReal solvetime; /* time at the conclusion of TSSolve() */
180: TSConvergedReason reason;
181: PetscBool errorifstepfailed;
182: PetscInt reject,max_reject;
183: TSExactFinalTimeOption exact_final_time;
185: PetscReal atol,rtol; /* Relative and absolute tolerance for local truncation error */
186: Vec vatol,vrtol; /* Relative and absolute tolerance in vector form */
187: PetscReal cfltime,cfltime_local;
189: /* ------------------- Default work-area management ------------------ */
190: PetscInt nwork;
191: Vec *work;
192: };
194: struct _TSAdaptOps {
195: PetscErrorCode (*choose)(TSAdapt,TS,PetscReal,PetscInt*,PetscReal*,PetscBool*,PetscReal*);
196: PetscErrorCode (*destroy)(TSAdapt);
197: PetscErrorCode (*reset)(TSAdapt);
198: PetscErrorCode (*view)(TSAdapt,PetscViewer);
199: PetscErrorCode (*setfromoptions)(PetscOptionItems*,TSAdapt);
200: PetscErrorCode (*load)(TSAdapt,PetscViewer);
201: };
203: struct _p_TSAdapt {
204: PETSCHEADER(struct _TSAdaptOps);
205: void *data;
206: PetscErrorCode (*checkstage)(TSAdapt,TS,PetscReal,Vec,PetscBool*);
207: struct {
208: PetscInt n; /* number of candidate schemes, including the one currently in use */
209: PetscBool inuse_set; /* the current scheme has been set */
210: const char *name[16]; /* name of the scheme */
211: PetscInt order[16]; /* classical order of each scheme */
212: PetscInt stageorder[16]; /* stage order of each scheme */
213: PetscReal ccfl[16]; /* stability limit relative to explicit Euler */
214: PetscReal cost[16]; /* relative measure of the amount of work required for each scheme */
215: } candidates;
216: PetscReal dt_min,dt_max;
217: PetscReal scale_solve_failed; /* Scale step by this factor if solver (linear or nonlinear) fails. */
218: PetscViewer monitor;
219: NormType wnormtype;
220: };
222: typedef struct _p_DMTS *DMTS;
223: typedef struct _DMTSOps *DMTSOps;
224: struct _DMTSOps {
225: TSRHSFunction rhsfunction;
226: TSRHSJacobian rhsjacobian;
228: TSIFunction ifunction;
229: PetscErrorCode (*ifunctionview)(void*,PetscViewer);
230: PetscErrorCode (*ifunctionload)(void**,PetscViewer);
232: TSIJacobian ijacobian;
233: PetscErrorCode (*ijacobianview)(void*,PetscViewer);
234: PetscErrorCode (*ijacobianload)(void**,PetscViewer);
236: TSI2Function i2function;
237: TSI2Jacobian i2jacobian;
239: TSSolutionFunction solution;
240: TSForcingFunction forcing;
242: PetscErrorCode (*destroy)(DMTS);
243: PetscErrorCode (*duplicate)(DMTS,DMTS);
244: };
246: struct _p_DMTS {
247: PETSCHEADER(struct _DMTSOps);
248: void *rhsfunctionctx;
249: void *rhsjacobianctx;
251: void *ifunctionctx;
252: void *ijacobianctx;
254: void *i2functionctx;
255: void *i2jacobianctx;
257: void *solutionctx;
258: void *forcingctx;
260: void *data;
262: /* This is NOT reference counted. The DM on which this context was first created is cached here to implement one-way
263: * copy-on-write. When DMGetDMTSWrite() sees a request using a different DM, it makes a copy. Thus, if a user
264: * only interacts directly with one level, e.g., using TSSetIFunction(), then coarse levels of a multilevel item
265: * integrator are built, then the user changes the routine with another call to TSSetIFunction(), it automatically
266: * propagates to all the levels. If instead, they get out a specific level and set the function on that level,
267: * subsequent changes to the original level will no longer propagate to that level.
268: */
269: DM originaldm;
270: };
272: PETSC_EXTERN PetscErrorCode DMGetDMTS(DM,DMTS*);
273: PETSC_EXTERN PetscErrorCode DMGetDMTSWrite(DM,DMTS*);
274: PETSC_EXTERN PetscErrorCode DMCopyDMTS(DM,DM);
275: PETSC_EXTERN PetscErrorCode DMTSView(DMTS,PetscViewer);
276: PETSC_EXTERN PetscErrorCode DMTSLoad(DMTS,PetscViewer);
277: PETSC_EXTERN PetscErrorCode DMTSCopy(DMTS,DMTS);
279: typedef enum {TSEVENT_NONE,TSEVENT_LOCATED_INTERVAL,TSEVENT_PROCESSING,TSEVENT_ZERO,TSEVENT_RESET_NEXTSTEP} TSEventStatus;
281: struct _n_TSEvent {
282: PetscScalar *fvalue; /* value of event function at the end of the step*/
283: PetscScalar *fvalue_prev; /* value of event function at start of the step (left end-point of event interval) */
284: PetscReal ptime_prev; /* time at step start (left end-point of event interval) */
285: 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) */
286: PetscReal ptime_right; /* time on the right end-point of the event interval */
287: PetscScalar *fvalue_right; /* value of event function at the right end-point of the event interval */
288: PetscInt *side; /* Used for detecting repetition of end-point, -1 => left, +1 => right */
289: PetscReal timestep_prev; /* previous time step */
290: PetscReal timestep_orig; /* initial time step */
291: PetscBool *zerocrossing; /* Flag to signal zero crossing detection */
292: PetscErrorCode (*eventhandler)(TS,PetscReal,Vec,PetscScalar*,void*); /* User event handler function */
293: PetscErrorCode (*postevent)(TS,PetscInt,PetscInt[],PetscReal,Vec,PetscBool,void*); /* User post event function */
294: void *ctx; /* User context for event handler and post even functions */
295: PetscInt *direction; /* Zero crossing direction: 1 -> Going positive, -1 -> Going negative, 0 -> Any */
296: PetscBool *terminate; /* 1 -> Terminate time stepping, 0 -> continue */
297: PetscInt nevents; /* Number of events to handle */
298: PetscInt nevents_zero; /* Number of event zero detected */
299: PetscInt *events_zero; /* List of events that have reached zero */
300: PetscReal *vtol; /* Vector tolerances for event zero check */
301: TSEventStatus status; /* Event status */
302: PetscInt iterctr; /* Iteration counter */
303: PetscViewer monitor;
304: /* Struct to record the events */
305: struct {
306: PetscInt ctr; /* recorder counter */
307: PetscReal *time; /* Event times */
308: PetscInt *stepnum; /* Step numbers */
309: PetscInt *nevents; /* Number of events occuring at the event times */
310: PetscInt **eventidx; /* Local indices of the events in the event list */
311: } recorder;
312: PetscInt recsize; /* Size of recorder stack */
313: };
315: PETSC_EXTERN PetscErrorCode TSEventInitialize(TSEvent,TS,PetscReal,Vec);
316: PETSC_EXTERN PetscErrorCode TSEventDestroy(TSEvent*);
317: PETSC_EXTERN PetscErrorCode TSEventHandler(TS);
318: PETSC_EXTERN PetscErrorCode TSAdjointEventHandler(TS);
320: PETSC_EXTERN PetscLogEvent TS_AdjointStep, TS_Step, TS_PseudoComputeTimeStep, TS_FunctionEval, TS_JacobianEval;
322: typedef enum {TS_STEP_INCOMPLETE, /* vec_sol, ptime, etc point to beginning of step */
323: TS_STEP_PENDING, /* vec_sol advanced, but step has not been accepted yet */
324: TS_STEP_COMPLETE /* step accepted and ptime, steps, etc have been advanced */
325: } TSStepStatus;
327: struct _n_TSMonitorLGCtx {
328: PetscDrawLG lg;
329: PetscInt howoften; /* when > 0 uses step % howoften, when negative only final solution plotted */
330: PetscInt ksp_its,snes_its;
331: char **names;
332: char **displaynames;
333: PetscInt ndisplayvariables;
334: PetscInt *displayvariables;
335: PetscReal *displayvalues;
336: PetscErrorCode (*transform)(void*,Vec,Vec*);
337: PetscErrorCode (*transformdestroy)(void*);
338: void *transformctx;
339: };
341: struct _n_TSMonitorEnvelopeCtx {
342: Vec max,min;
343: };
345: PETSC_EXTERN PetscLogEvent TSTrajectory_Set, TSTrajectory_Get, TSTrajectory_DiskWrite, TSTrajectory_DiskRead;
347: #endif