Actual source code: tsimpl.h
petsc-3.8.4 2018-03-24
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 TSGLLERegisterAll(void);
28: PETSC_EXTERN PetscErrorCode TSGLLEAdaptRegisterAll(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 (*forwardsetup)(TS);
53: PetscErrorCode (*forwardstep)(TS);
54: PetscErrorCode (*forwardintegral)(TS);
55: PetscErrorCode (*getsolutioncomponents)(TS,PetscInt*,Vec*);
56: PetscErrorCode (*getauxsolution)(TS,Vec*);
57: PetscErrorCode (*gettimeerror)(TS,PetscInt,Vec*);
58: PetscErrorCode (*settimeerror)(TS,Vec);
59: PetscErrorCode (*startingmethod) (TS);
60: };
62: /*
63: TSEvent - Abstract object to handle event monitoring
64: */
65: typedef struct _n_TSEvent *TSEvent;
67: typedef struct _TSTrajectoryOps *TSTrajectoryOps;
69: struct _TSTrajectoryOps {
70: PetscErrorCode (*view)(TSTrajectory,PetscViewer);
71: PetscErrorCode (*destroy)(TSTrajectory);
72: PetscErrorCode (*set)(TSTrajectory,TS,PetscInt,PetscReal,Vec);
73: PetscErrorCode (*get)(TSTrajectory,TS,PetscInt,PetscReal*);
74: PetscErrorCode (*setfromoptions)(PetscOptionItems*,TSTrajectory);
75: PetscErrorCode (*setup)(TSTrajectory,TS);
76: };
78: struct _p_TSTrajectory {
79: PETSCHEADER(struct _TSTrajectoryOps);
80: PetscViewer monitor;
81: PetscInt setupcalled; /* true if setup has been called */
82: PetscInt recomps; /* counter for recomputations in the adjoint run */
83: PetscInt diskreads,diskwrites; /* counters for disk checkpoint reads and writes */
84: char **names; /* the name of each variable; each process has only the local names */
85: PetscErrorCode (*transform)(void*,Vec,Vec*);
86: PetscErrorCode (*transformdestroy)(void*);
87: void* transformctx;
88: void *data;
89: };
91: struct _p_TS {
92: PETSCHEADER(struct _TSOps);
93: TSProblemType problem_type;
94: TSEquationType equation_type;
96: DM dm;
97: Vec vec_sol; /* solution vector in first and second order equations */
98: Vec vec_dot; /* time derivative vector in second order equations */
99: TSAdapt adapt;
100: TSAdaptType default_adapt_type;
101: TSEvent event;
103: /* ---------------- User (or PETSc) Provided stuff ---------------------*/
104: PetscErrorCode (*monitor[MAXTSMONITORS])(TS,PetscInt,PetscReal,Vec,void*);
105: PetscErrorCode (*monitordestroy[MAXTSMONITORS])(void**);
106: void *monitorcontext[MAXTSMONITORS];
107: PetscInt numbermonitors;
108: PetscErrorCode (*adjointmonitor[MAXTSMONITORS])(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*);
109: PetscErrorCode (*adjointmonitordestroy[MAXTSMONITORS])(void**);
110: void *adjointmonitorcontext[MAXTSMONITORS];
111: PetscInt numberadjointmonitors;
113: PetscErrorCode (*prestep)(TS);
114: PetscErrorCode (*prestage)(TS,PetscReal);
115: PetscErrorCode (*poststage)(TS,PetscReal,PetscInt,Vec*);
116: PetscErrorCode (*postevaluate)(TS);
117: PetscErrorCode (*poststep)(TS);
118: PetscErrorCode (*functiondomainerror)(TS,PetscReal,Vec,PetscBool*);
120: /* ---------------------- Sensitivity Analysis support -----------------*/
121: TSTrajectory trajectory; /* All solutions are kept here for the entire time integration process */
122: Vec *vecs_sensi; /* one vector for each cost function */
123: Vec *vecs_sensip;
124: PetscInt numcost; /* number of cost functions */
125: Vec vec_costintegral;
126: PetscInt adjointsetupcalled;
127: PetscInt adjoint_steps;
128: PetscInt adjoint_max_steps;
129: PetscBool adjoint_solve; /* immediately call TSAdjointSolve() after TSSolve() is complete */
130: PetscBool costintegralfwd; /* cost integral is evaluated in the forward run if true */
131: Vec vec_costintegrand; /* workspace for Adjoint computations */
132: Mat Jacp;
133: void *rhsjacobianpctx;
134: void *costintegrandctx;
135: Vec *vecs_drdy;
136: Vec *vecs_drdp;
138: PetscErrorCode (*rhsjacobianp)(TS,PetscReal,Vec,Mat,void*);
139: PetscErrorCode (*costintegrand)(TS,PetscReal,Vec,Vec,void*);
140: PetscErrorCode (*drdyfunction)(TS,PetscReal,Vec,Vec*,void*);
141: PetscErrorCode (*drdpfunction)(TS,PetscReal,Vec,Vec*,void*);
143: /* specific to forward sensitivity analysis */
144: Vec *vecs_fwdsensipacked; /* packed vector array for forward sensitivitis */
145: Vec *vecs_integral_sensi; /* one vector for each integral */
146: Vec *vecs_integral_sensip; /* one vector for each integral */
147: PetscInt num_parameters;
148: PetscInt num_initialvalues;
149: Vec *vecs_jacp;
150: void *vecsrhsjacobianpctx;
151: PetscInt forwardsetupcalled;
152: PetscBool forward_solve;
153: PetscErrorCode (*vecsrhsjacobianp)(TS,PetscReal,Vec,Vec*,void*);
155: /* ---------------------- IMEX support ---------------------------------*/
156: /* These extra slots are only used when the user provides both Implicit and RHS */
157: Mat Arhs; /* Right hand side matrix */
158: Mat Brhs; /* Right hand side preconditioning matrix */
159: Vec Frhs; /* Right hand side function value */
161: /* This is a general caching scheme to avoid recomputing the Jacobian at a place that has been previously been evaluated.
162: * The present use case is that TSComputeRHSFunctionLinear() evaluates the Jacobian once and we don't want it to be immeditely re-evaluated.
163: */
164: struct {
165: PetscReal time; /* The time at which the matrices were last evaluated */
166: PetscObjectId Xid; /* Unique ID of solution vector at which the Jacobian was last evaluated */
167: PetscObjectState Xstate; /* State of the solution vector */
168: MatStructure mstructure; /* The structure returned */
169: /* Flag to unshift Jacobian before calling the IJacobian or RHSJacobian functions. This is useful
170: * if the user would like to reuse (part of) the Jacobian from the last evaluation. */
171: PetscBool reuse;
172: PetscReal scale,shift;
173: } rhsjacobian;
175: struct {
176: PetscReal shift; /* The derivative of the lhs wrt to Xdot */
177: } ijacobian;
179: /* --------------------Nonlinear Iteration------------------------------*/
180: SNES snes;
181: PetscBool usessnes; /* Flag set by each TSType to indicate if the type actually uses a SNES;
182: this works around the design flaw that a SNES is ALWAYS created with TS even when it is not needed.*/
183: PetscInt ksp_its; /* total number of linear solver iterations */
184: PetscInt snes_its; /* total number of nonlinear solver iterations */
185: PetscInt num_snes_failures;
186: PetscInt max_snes_failures;
188: /* --- Data that is unique to each particular solver --- */
189: PetscInt setupcalled; /* true if setup has been called */
190: void *data; /* implementationspecific data */
191: void *user; /* user context */
193: /* ------------------ Parameters -------------------------------------- */
194: PetscInt max_steps; /* max number of steps */
195: PetscReal max_time; /* max time allowed */
197: /* --------------------------------------------------------------------- */
199: PetscBool steprollback; /* flag to indicate that the step was rolled back */
200: PetscBool steprestart; /* flag to indicate that the timestepper has to discard any history and restart */
201: PetscInt steps; /* steps taken so far in all successive calls to TSSolve() */
202: PetscReal ptime; /* time at the start of the current step (stage time is internal if it exists) */
203: PetscReal time_step; /* current time increment */
204: PetscReal ptime_prev; /* time at the start of the previous step */
205: PetscReal ptime_prev_rollback; /* time at the start of the 2nd previous step to recover from rollback */
206: PetscReal solvetime; /* time at the conclusion of TSSolve() */
208: TSConvergedReason reason;
209: PetscBool errorifstepfailed;
210: PetscInt reject,max_reject;
211: TSExactFinalTimeOption exact_final_time;
213: PetscReal atol,rtol; /* Relative and absolute tolerance for local truncation error */
214: Vec vatol,vrtol; /* Relative and absolute tolerance in vector form */
215: PetscReal cfltime,cfltime_local;
217: /* ------------------- Default work-area management ------------------ */
218: PetscInt nwork;
219: Vec *work;
220: };
222: struct _TSAdaptOps {
223: PetscErrorCode (*choose)(TSAdapt,TS,PetscReal,PetscInt*,PetscReal*,PetscBool*,PetscReal*,PetscReal*,PetscReal*);
224: PetscErrorCode (*destroy)(TSAdapt);
225: PetscErrorCode (*reset)(TSAdapt);
226: PetscErrorCode (*view)(TSAdapt,PetscViewer);
227: PetscErrorCode (*setfromoptions)(PetscOptionItems*,TSAdapt);
228: PetscErrorCode (*load)(TSAdapt,PetscViewer);
229: };
231: struct _p_TSAdapt {
232: PETSCHEADER(struct _TSAdaptOps);
233: void *data;
234: PetscErrorCode (*checkstage)(TSAdapt,TS,PetscReal,Vec,PetscBool*);
235: struct {
236: PetscInt n; /* number of candidate schemes, including the one currently in use */
237: PetscBool inuse_set; /* the current scheme has been set */
238: const char *name[16]; /* name of the scheme */
239: PetscInt order[16]; /* classical order of each scheme */
240: PetscInt stageorder[16]; /* stage order of each scheme */
241: PetscReal ccfl[16]; /* stability limit relative to explicit Euler */
242: PetscReal cost[16]; /* relative measure of the amount of work required for each scheme */
243: } candidates;
244: PetscBool always_accept;
245: PetscReal safety; /* safety factor relative to target error/stability goal */
246: PetscReal reject_safety; /* extra safety factor if the last step was rejected */
247: PetscReal clip[2]; /* admissible time step decrease/increase factors */
248: PetscReal dt_min,dt_max; /* admissible minimum and maximum time step */
249: PetscReal scale_solve_failed; /* scale step by this factor if solver (linear or nonlinear) fails. */
250: NormType wnormtype;
251: PetscViewer monitor;
252: PetscInt timestepjustincreased;
253: };
255: typedef struct _p_DMTS *DMTS;
256: typedef struct _DMTSOps *DMTSOps;
257: struct _DMTSOps {
258: TSRHSFunction rhsfunction;
259: TSRHSJacobian rhsjacobian;
261: TSIFunction ifunction;
262: PetscErrorCode (*ifunctionview)(void*,PetscViewer);
263: PetscErrorCode (*ifunctionload)(void**,PetscViewer);
265: TSIJacobian ijacobian;
266: PetscErrorCode (*ijacobianview)(void*,PetscViewer);
267: PetscErrorCode (*ijacobianload)(void**,PetscViewer);
269: TSI2Function i2function;
270: TSI2Jacobian i2jacobian;
272: TSSolutionFunction solution;
273: TSForcingFunction forcing;
275: PetscErrorCode (*destroy)(DMTS);
276: PetscErrorCode (*duplicate)(DMTS,DMTS);
277: };
279: struct _p_DMTS {
280: PETSCHEADER(struct _DMTSOps);
281: void *rhsfunctionctx;
282: void *rhsjacobianctx;
284: void *ifunctionctx;
285: void *ijacobianctx;
287: void *i2functionctx;
288: void *i2jacobianctx;
290: void *solutionctx;
291: void *forcingctx;
293: void *data;
295: /* This is NOT reference counted. The DM on which this context was first created is cached here to implement one-way
296: * copy-on-write. When DMGetDMTSWrite() sees a request using a different DM, it makes a copy. Thus, if a user
297: * only interacts directly with one level, e.g., using TSSetIFunction(), then coarse levels of a multilevel item
298: * integrator are built, then the user changes the routine with another call to TSSetIFunction(), it automatically
299: * propagates to all the levels. If instead, they get out a specific level and set the function on that level,
300: * subsequent changes to the original level will no longer propagate to that level.
301: */
302: DM originaldm;
303: };
305: PETSC_EXTERN PetscErrorCode DMGetDMTS(DM,DMTS*);
306: PETSC_EXTERN PetscErrorCode DMGetDMTSWrite(DM,DMTS*);
307: PETSC_EXTERN PetscErrorCode DMCopyDMTS(DM,DM);
308: PETSC_EXTERN PetscErrorCode DMTSView(DMTS,PetscViewer);
309: PETSC_EXTERN PetscErrorCode DMTSLoad(DMTS,PetscViewer);
310: PETSC_EXTERN PetscErrorCode DMTSCopy(DMTS,DMTS);
312: typedef enum {TSEVENT_NONE,TSEVENT_LOCATED_INTERVAL,TSEVENT_PROCESSING,TSEVENT_ZERO,TSEVENT_RESET_NEXTSTEP} TSEventStatus;
314: struct _n_TSEvent {
315: PetscScalar *fvalue; /* value of event function at the end of the step*/
316: PetscScalar *fvalue_prev; /* value of event function at start of the step (left end-point of event interval) */
317: PetscReal ptime_prev; /* time at step start (left end-point of event interval) */
318: 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) */
319: PetscReal ptime_right; /* time on the right end-point of the event interval */
320: PetscScalar *fvalue_right; /* value of event function at the right end-point of the event interval */
321: PetscInt *side; /* Used for detecting repetition of end-point, -1 => left, +1 => right */
322: PetscReal timestep_prev; /* previous time step */
323: PetscReal timestep_orig; /* initial time step */
324: PetscBool *zerocrossing; /* Flag to signal zero crossing detection */
325: PetscErrorCode (*eventhandler)(TS,PetscReal,Vec,PetscScalar*,void*); /* User event handler function */
326: PetscErrorCode (*postevent)(TS,PetscInt,PetscInt[],PetscReal,Vec,PetscBool,void*); /* User post event function */
327: void *ctx; /* User context for event handler and post even functions */
328: PetscInt *direction; /* Zero crossing direction: 1 -> Going positive, -1 -> Going negative, 0 -> Any */
329: PetscBool *terminate; /* 1 -> Terminate time stepping, 0 -> continue */
330: PetscInt nevents; /* Number of events to handle */
331: PetscInt nevents_zero; /* Number of event zero detected */
332: PetscInt *events_zero; /* List of events that have reached zero */
333: PetscReal *vtol; /* Vector tolerances for event zero check */
334: TSEventStatus status; /* Event status */
335: PetscInt iterctr; /* Iteration counter */
336: PetscViewer monitor;
337: /* Struct to record the events */
338: struct {
339: PetscInt ctr; /* recorder counter */
340: PetscReal *time; /* Event times */
341: PetscInt *stepnum; /* Step numbers */
342: PetscInt *nevents; /* Number of events occuring at the event times */
343: PetscInt **eventidx; /* Local indices of the events in the event list */
344: } recorder;
345: PetscInt recsize; /* Size of recorder stack */
346: PetscInt refct; /* reference count */
347: };
349: PETSC_EXTERN PetscErrorCode TSEventInitialize(TSEvent,TS,PetscReal,Vec);
350: PETSC_EXTERN PetscErrorCode TSEventDestroy(TSEvent*);
351: PETSC_EXTERN PetscErrorCode TSEventHandler(TS);
352: PETSC_EXTERN PetscErrorCode TSAdjointEventHandler(TS);
354: PETSC_EXTERN PetscLogEvent TS_AdjointStep, TS_Step, TS_PseudoComputeTimeStep, TS_FunctionEval, TS_JacobianEval, TS_ForwardStep;
356: typedef enum {TS_STEP_INCOMPLETE, /* vec_sol, ptime, etc point to beginning of step */
357: TS_STEP_PENDING, /* vec_sol advanced, but step has not been accepted yet */
358: TS_STEP_COMPLETE /* step accepted and ptime, steps, etc have been advanced */
359: } TSStepStatus;
361: struct _n_TSMonitorLGCtx {
362: PetscDrawLG lg;
363: PetscBool semilogy;
364: PetscInt howoften; /* when > 0 uses step % howoften, when negative only final solution plotted */
365: PetscInt ksp_its,snes_its;
366: char **names;
367: char **displaynames;
368: PetscInt ndisplayvariables;
369: PetscInt *displayvariables;
370: PetscReal *displayvalues;
371: PetscErrorCode (*transform)(void*,Vec,Vec*);
372: PetscErrorCode (*transformdestroy)(void*);
373: void *transformctx;
374: };
376: struct _n_TSMonitorEnvelopeCtx {
377: Vec max,min;
378: };
380: /*
381: Checks if the user provide a TSSetIFunction() but an explicit method is called; generate an error in that case
382: */
383: PETSC_STATIC_INLINE PetscErrorCode TSCheckImplicitTerm(TS ts)
384: {
385: TSIFunction ifunction;
386: DM dm;
387: PetscErrorCode ierr;
390: TSGetDM(ts,&dm);
391: DMTSGetIFunction(dm,&ifunction,NULL);
392: 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()");
393: return(0);
394: }
396: PETSC_EXTERN PetscLogEvent TSTrajectory_Set, TSTrajectory_Get, TSTrajectory_DiskWrite, TSTrajectory_DiskRead;
398: #endif