Actual source code: petsc-tsimpl.h
petsc-dev 2014-02-02
2: #ifndef __TSIMPL_H
5: #include <petscts.h>
6: #include <petsc-private/petscimpl.h>
8: /*
9: Timesteping context.
10: General DAE: F(t,U,U_t) = 0, required Jacobian is G'(U) where G(U) = F(t,U,U0+a*U)
11: General ODE: U_t = F(t,U) <-- the right-hand-side function
12: Linear ODE: U_t = A(t) U <-- the right-hand-side matrix
13: Linear (no time) ODE: U_t = A U <-- the right-hand-side matrix
14: */
16: /*
17: Maximum number of monitors you can run with a single TS
18: */
19: #define MAXTSMONITORS 5
21: typedef struct _TSOps *TSOps;
23: struct _TSOps {
24: PetscErrorCode (*snesfunction)(SNES,Vec,Vec,TS);
25: PetscErrorCode (*snesjacobian)(SNES,Vec,Mat*,Mat*,MatStructure*,TS);
26: PetscErrorCode (*setup)(TS);
27: PetscErrorCode (*step)(TS);
28: PetscErrorCode (*solve)(TS);
29: PetscErrorCode (*interpolate)(TS,PetscReal,Vec);
30: PetscErrorCode (*evaluatestep)(TS,PetscInt,Vec,PetscBool*);
31: PetscErrorCode (*setfromoptions)(TS);
32: PetscErrorCode (*destroy)(TS);
33: PetscErrorCode (*view)(TS,PetscViewer);
34: PetscErrorCode (*reset)(TS);
35: PetscErrorCode (*linearstability)(TS,PetscReal,PetscReal,PetscReal*,PetscReal*);
36: PetscErrorCode (*load)(TS,PetscViewer);
37: };
39: struct _p_TS {
40: PETSCHEADER(struct _TSOps);
41: DM dm;
42: TSProblemType problem_type;
43: Vec vec_sol;
44: TSAdapt adapt;
46: /* ---------------- User (or PETSc) Provided stuff ---------------------*/
47: PetscErrorCode (*monitor[MAXTSMONITORS])(TS,PetscInt,PetscReal,Vec,void*); /* returns control to user after */
48: PetscErrorCode (*monitordestroy[MAXTSMONITORS])(void**);
49: void *monitorcontext[MAXTSMONITORS]; /* residual calculation, allows user */
50: PetscInt numbermonitors; /* to, for instance, print residual norm, etc. */
52: PetscErrorCode (*prestep)(TS);
53: PetscErrorCode (*prestage)(TS,PetscReal);
54: PetscErrorCode (*poststage)(TS,PetscReal,PetscInt,Vec*);
55: PetscErrorCode (*poststep)(TS);
57: /* ---------------------- IMEX support ---------------------------------*/
58: /* These extra slots are only used when the user provides both Implicit and RHS */
59: Mat Arhs; /* Right hand side matrix */
60: Mat Brhs; /* Right hand side preconditioning matrix */
61: Vec Frhs; /* Right hand side function value */
63: /* This is a general caching scheme to avoid recomputing the Jacobian at a place that has been previously been evaluated.
64: * The present use case is that TSComputeRHSFunctionLinear() evaluates the Jacobian once and we don't want it to be immeditely re-evaluated.
65: */
66: struct {
67: PetscReal time; /* The time at which the matrices were last evaluated */
68: Vec X; /* Solution vector at which the Jacobian was last evaluated */
69: PetscObjectState Xstate; /* State of the solution vector */
70: MatStructure mstructure; /* The structure returned */
71: /* Flag to unshift Jacobian before calling the IJacobian or RHSJacobian functions. This is useful
72: * if the user would like to reuse (part of) the Jacobian from the last evaluation. */
73: PetscBool reuse;
74: PetscReal scale,shift;
75: } rhsjacobian;
77: struct {
78: PetscReal shift; /* The derivative of the lhs wrt to Xdot */
79: } ijacobian;
81: /* ---------------------Nonlinear Iteration------------------------------*/
82: SNES snes;
84: /* --- Data that is unique to each particular solver --- */
85: PetscInt setupcalled; /* true if setup has been called */
86: void *data; /* implementationspecific data */
87: void *user; /* user context */
89: /* ------------------ Parameters -------------------------------------- */
90: PetscInt max_steps; /* max number of steps */
91: PetscReal max_time; /* max time allowed */
92: PetscReal time_step; /* current/completed time increment */
93: PetscReal time_step_prev; /* previous time step */
95: /*
96: these are temporary to support increasing the time step if nonlinear solver convergence remains good
97: and time_step was previously cut due to failed nonlinear solver
98: */
99: PetscReal time_step_orig; /* original time step requested by user */
100: PetscInt time_steps_since_decrease; /* number of timesteps since timestep was decreased due to lack of convergence */
101: /* ----------------------------------------------------------------------------------------------------------------*/
103: PetscInt steps; /* steps taken so far */
104: PetscReal ptime; /* time at the start of the current step (stage time is internal if it exists) */
105: PetscReal solvetime; /* time at the conclusion of TSSolve() */
106: PetscInt ksp_its; /* total number of linear solver iterations */
107: PetscInt snes_its; /* total number of nonlinear solver iterations */
109: PetscInt num_snes_failures;
110: PetscInt max_snes_failures;
111: TSConvergedReason reason;
112: TSEquationType equation_type;
113: PetscBool errorifstepfailed;
114: TSExactFinalTimeOption exact_final_time;
115: PetscBool retain_stages;
116: PetscInt reject,max_reject;
118: PetscReal atol,rtol; /* Relative and absolute tolerance for local truncation error */
119: Vec vatol,vrtol; /* Relative and absolute tolerance in vector form */
120: PetscReal cfltime,cfltime_local;
122: /* ------------------- Default work-area management ------------------ */
123: PetscInt nwork;
124: Vec *work;
125: };
127: struct _TSAdaptOps {
128: PetscErrorCode (*choose)(TSAdapt,TS,PetscReal,PetscInt*,PetscReal*,PetscBool*,PetscReal*);
129: PetscErrorCode (*checkstage)(TSAdapt,TS,PetscBool*);
130: PetscErrorCode (*destroy)(TSAdapt);
131: PetscErrorCode (*view)(TSAdapt,PetscViewer);
132: PetscErrorCode (*setfromoptions)(TSAdapt);
133: PetscErrorCode (*load)(TSAdapt,PetscViewer);
134: };
136: struct _p_TSAdapt {
137: PETSCHEADER(struct _TSAdaptOps);
138: void *data;
139: struct {
140: PetscInt n; /* number of candidate schemes, including the one currently in use */
141: PetscBool inuse_set; /* the current scheme has been set */
142: const char *name[16]; /* name of the scheme */
143: PetscInt order[16]; /* classical order of each scheme */
144: PetscInt stageorder[16]; /* stage order of each scheme */
145: PetscReal ccfl[16]; /* stability limit relative to explicit Euler */
146: PetscReal cost[16]; /* relative measure of the amount of work required for each scheme */
147: } candidates;
148: PetscReal dt_min,dt_max;
149: PetscReal scale_solve_failed; /* Scale step by this factor if solver (linear or nonlinear) fails. */
150: PetscViewer monitor;
151: };
153: typedef struct _p_DMTS *DMTS;
154: typedef struct _DMTSOps *DMTSOps;
155: struct _DMTSOps {
156: TSRHSFunction rhsfunction;
157: TSRHSJacobian rhsjacobian;
159: TSIFunction ifunction;
160: PetscErrorCode (*ifunctionview)(void*,PetscViewer);
161: PetscErrorCode (*ifunctionload)(void**,PetscViewer);
163: TSIJacobian ijacobian;
164: PetscErrorCode (*ijacobianview)(void*,PetscViewer);
165: PetscErrorCode (*ijacobianload)(void**,PetscViewer);
167: TSSolutionFunction solution;
168: PetscErrorCode (*forcing)(TS,PetscReal,Vec,void*);
170: PetscErrorCode (*destroy)(DMTS);
171: PetscErrorCode (*duplicate)(DMTS,DMTS);
172: };
174: struct _p_DMTS {
175: PETSCHEADER(struct _DMTSOps);
176: void *rhsfunctionctx;
177: void *rhsjacobianctx;
179: void *ifunctionctx;
180: void *ijacobianctx;
182: void *solutionctx;
183: void *forcingctx;
185: void *data;
187: /* This is NOT reference counted. The DM on which this context was first created is cached here to implement one-way
188: * copy-on-write. When DMGetDMTSWrite() sees a request using a different DM, it makes a copy. Thus, if a user
189: * only interacts directly with one level, e.g., using TSSetIFunction(), then coarse levels of a multilevel item
190: * integrator are built, then the user changes the routine with another call to TSSetIFunction(), it automatically
191: * propagates to all the levels. If instead, they get out a specific level and set the function on that level,
192: * subsequent changes to the original level will no longer propagate to that level.
193: */
194: DM originaldm;
195: };
197: PETSC_EXTERN PetscErrorCode DMGetDMTS(DM,DMTS*);
198: PETSC_EXTERN PetscErrorCode DMGetDMTSWrite(DM,DMTS*);
199: PETSC_EXTERN PetscErrorCode DMCopyDMTS(DM,DM);
200: PETSC_EXTERN PetscErrorCode DMTSView(DMTS,PetscViewer);
201: PETSC_EXTERN PetscErrorCode DMTSLoad(DMTS,PetscViewer);
204: PETSC_EXTERN PetscLogEvent TS_Step, TS_PseudoComputeTimeStep, TS_FunctionEval, TS_JacobianEval;
206: typedef enum {TS_STEP_INCOMPLETE, /* vec_sol, ptime, etc point to beginning of step */
207: TS_STEP_PENDING, /* vec_sol advanced, but step has not been accepted yet */
208: TS_STEP_COMPLETE /* step accepted and ptime, steps, etc have been advanced */
209: } TSStepStatus;
211: struct _n_TSMonitorLGCtx {
212: PetscDrawLG lg;
213: PetscInt howoften; /* when > 0 uses step % howoften, when negative only final solution plotted */
214: PetscInt ksp_its,snes_its;
215: };
217: #endif