Actual source code: tsimpl.h

petsc-3.9.2 2018-05-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 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:   PetscBool      keepfiles;               /* keep the files generated during the run after the run is complete */
 86:   char           *dirname,*filetemplate;  /* directory name and file name template for disk checkpoints */
 87:   char           *dirfiletemplate;        /* complete directory and file name template for disk checkpoints */
 88:   PetscErrorCode (*transform)(void*,Vec,Vec*);
 89:   PetscErrorCode (*transformdestroy)(void*);
 90:   void*          transformctx;
 91:   void           *data;
 92: };

 94: struct _p_TS {
 95:   PETSCHEADER(struct _TSOps);
 96:   TSProblemType  problem_type;
 97:   TSEquationType equation_type;

 99:   DM             dm;
100:   Vec            vec_sol; /* solution vector in first and second order equations */
101:   Vec            vec_dot; /* time derivative vector in second order equations */
102:   TSAdapt        adapt;
103:   TSAdaptType    default_adapt_type;
104:   TSEvent        event;

106:   /* ---------------- User (or PETSc) Provided stuff ---------------------*/
107:   PetscErrorCode (*monitor[MAXTSMONITORS])(TS,PetscInt,PetscReal,Vec,void*);
108:   PetscErrorCode (*monitordestroy[MAXTSMONITORS])(void**);
109:   void *monitorcontext[MAXTSMONITORS];
110:   PetscInt  numbermonitors;
111:   PetscErrorCode (*adjointmonitor[MAXTSMONITORS])(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*);
112:   PetscErrorCode (*adjointmonitordestroy[MAXTSMONITORS])(void**);
113:   void *adjointmonitorcontext[MAXTSMONITORS];
114:   PetscInt  numberadjointmonitors;

116:   PetscErrorCode (*prestep)(TS);
117:   PetscErrorCode (*prestage)(TS,PetscReal);
118:   PetscErrorCode (*poststage)(TS,PetscReal,PetscInt,Vec*);
119:   PetscErrorCode (*postevaluate)(TS);
120:   PetscErrorCode (*poststep)(TS);
121:   PetscErrorCode (*functiondomainerror)(TS,PetscReal,Vec,PetscBool*);

123:   /* ---------------------- Sensitivity Analysis support -----------------*/
124:   TSTrajectory trajectory;          /* All solutions are kept here for the entire time integration process */
125:   Vec       *vecs_sensi;            /* one vector for each cost function */
126:   Vec       *vecs_sensip;
127:   PetscInt  numcost;                /* number of cost functions */
128:   Vec       vec_costintegral;
129:   PetscInt  adjointsetupcalled;
130:   PetscInt  adjoint_steps;
131:   PetscInt  adjoint_max_steps;
132:   PetscBool adjoint_solve;          /* immediately call TSAdjointSolve() after TSSolve() is complete */
133:   PetscBool costintegralfwd;        /* cost integral is evaluated in the forward run if true */
134:   Vec       vec_costintegrand;      /* workspace for Adjoint computations */
135:   Mat       Jacp;
136:   void      *rhsjacobianpctx;
137:   void      *costintegrandctx;
138:   Vec       *vecs_drdy;
139:   Vec       *vecs_drdp;

141:   PetscErrorCode (*rhsjacobianp)(TS,PetscReal,Vec,Mat,void*);
142:   PetscErrorCode (*costintegrand)(TS,PetscReal,Vec,Vec,void*);
143:   PetscErrorCode (*drdyfunction)(TS,PetscReal,Vec,Vec*,void*);
144:   PetscErrorCode (*drdpfunction)(TS,PetscReal,Vec,Vec*,void*);

146:   /* specific to forward sensitivity analysis */
147:   Mat       mat_sensip;              /* matrix storing forward sensitivities */
148:   Vec       *vecs_integral_sensip;   /* one vector for each integral */
149:   PetscInt  num_parameters;
150:   PetscInt  num_initialvalues;
151:   void      *vecsrhsjacobianpctx;
152:   PetscInt  forwardsetupcalled;
153:   PetscBool forward_solve;
154:   PetscErrorCode (*vecsrhsjacobianp)(TS,PetscReal,Vec,Vec*,void*);

156:   /* ---------------------- IMEX support ---------------------------------*/
157:   /* These extra slots are only used when the user provides both Implicit and RHS */
158:   Mat Arhs;     /* Right hand side matrix */
159:   Mat Brhs;     /* Right hand side preconditioning matrix */
160:   Vec Frhs;     /* Right hand side function value */

162:   /* This is a general caching scheme to avoid recomputing the Jacobian at a place that has been previously been evaluated.
163:    * The present use case is that TSComputeRHSFunctionLinear() evaluates the Jacobian once and we don't want it to be immeditely re-evaluated.
164:    */
165:   struct {
166:     PetscReal        time;          /* The time at which the matrices were last evaluated */
167:     PetscObjectId    Xid;           /* Unique ID of solution vector at which the Jacobian was last evaluated */
168:     PetscObjectState Xstate;        /* State of the solution vector */
169:     MatStructure     mstructure;    /* The structure returned */
170:     /* Flag to unshift Jacobian before calling the IJacobian or RHSJacobian functions.  This is useful
171:      * if the user would like to reuse (part of) the Jacobian from the last evaluation. */
172:     PetscBool        reuse;
173:     PetscReal        scale,shift;
174:   } rhsjacobian;

176:   struct {
177:     PetscReal shift;            /* The derivative of the lhs wrt to Xdot */
178:   } ijacobian;

180:   /* --------------------Nonlinear Iteration------------------------------*/
181:   SNES     snes;
182:   PetscBool usessnes;   /* Flag set by each TSType to indicate if the type actually uses a SNES;
183:                            this works around the design flaw that a SNES is ALWAYS created with TS even when it is not needed.*/
184:   PetscInt ksp_its;                /* total number of linear solver iterations */
185:   PetscInt snes_its;               /* total number of nonlinear solver iterations */
186:   PetscInt num_snes_failures;
187:   PetscInt max_snes_failures;

189:   /* --- Data that is unique to each particular solver --- */
190:   PetscInt setupcalled;             /* true if setup has been called */
191:   void     *data;                   /* implementationspecific data */
192:   void     *user;                   /* user context */

194:   /* ------------------  Parameters -------------------------------------- */
195:   PetscInt  max_steps;              /* max number of steps */
196:   PetscReal max_time;               /* max time allowed */

198:   /* --------------------------------------------------------------------- */

200:   PetscBool steprollback;           /* flag to indicate that the step was rolled back */
201:   PetscBool steprestart;            /* flag to indicate that the timestepper has to discard any history and restart */
202:   PetscInt  steps;                  /* steps taken so far in all successive calls to TSSolve() */
203:   PetscReal ptime;                  /* time at the start of the current step (stage time is internal if it exists) */
204:   PetscReal time_step;              /* current time increment */
205:   PetscReal ptime_prev;             /* time at the start of the previous step */
206:   PetscReal ptime_prev_rollback;    /* time at the start of the 2nd previous step to recover from rollback */
207:   PetscReal solvetime;              /* time at the conclusion of TSSolve() */

209:   TSConvergedReason reason;
210:   PetscBool errorifstepfailed;
211:   PetscInt  reject,max_reject;
212:   TSExactFinalTimeOption exact_final_time;

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

218:   PetscBool testjacobian;
219:   PetscBool testjacobiantranspose;
220:   /* ------------------- Default work-area management ------------------ */
221:   PetscInt nwork;
222:   Vec      *work;
223: };

225: struct _TSAdaptOps {
226:   PetscErrorCode (*choose)(TSAdapt,TS,PetscReal,PetscInt*,PetscReal*,PetscBool*,PetscReal*,PetscReal*,PetscReal*);
227:   PetscErrorCode (*destroy)(TSAdapt);
228:   PetscErrorCode (*reset)(TSAdapt);
229:   PetscErrorCode (*view)(TSAdapt,PetscViewer);
230:   PetscErrorCode (*setfromoptions)(PetscOptionItems*,TSAdapt);
231:   PetscErrorCode (*load)(TSAdapt,PetscViewer);
232: };

234: struct _p_TSAdapt {
235:   PETSCHEADER(struct _TSAdaptOps);
236:   void *data;
237:   PetscErrorCode (*checkstage)(TSAdapt,TS,PetscReal,Vec,PetscBool*);
238:   struct {
239:     PetscInt   n;                /* number of candidate schemes, including the one currently in use */
240:     PetscBool  inuse_set;        /* the current scheme has been set */
241:     const char *name[16];        /* name of the scheme */
242:     PetscInt   order[16];        /* classical order of each scheme */
243:     PetscInt   stageorder[16];   /* stage order of each scheme */
244:     PetscReal  ccfl[16];         /* stability limit relative to explicit Euler */
245:     PetscReal  cost[16];         /* relative measure of the amount of work required for each scheme */
246:   } candidates;
247:   PetscBool   always_accept;
248:   PetscReal   safety;             /* safety factor relative to target error/stability goal */
249:   PetscReal   reject_safety;      /* extra safety factor if the last step was rejected */
250:   PetscReal   clip[2];            /* admissible time step decrease/increase factors */
251:   PetscReal   dt_min,dt_max;      /* admissible minimum and maximum time step */
252:   PetscReal   scale_solve_failed; /* scale step by this factor if solver (linear or nonlinear) fails. */
253:   NormType    wnormtype;
254:   PetscViewer monitor;
255:   PetscInt    timestepjustincreased;
256: };

258: typedef struct _p_DMTS *DMTS;
259: typedef struct _DMTSOps *DMTSOps;
260: struct _DMTSOps {
261:   TSRHSFunction rhsfunction;
262:   TSRHSJacobian rhsjacobian;

264:   TSIFunction ifunction;
265:   PetscErrorCode (*ifunctionview)(void*,PetscViewer);
266:   PetscErrorCode (*ifunctionload)(void**,PetscViewer);

268:   TSIJacobian ijacobian;
269:   PetscErrorCode (*ijacobianview)(void*,PetscViewer);
270:   PetscErrorCode (*ijacobianload)(void**,PetscViewer);

272:   TSI2Function i2function;
273:   TSI2Jacobian i2jacobian;

275:   TSSolutionFunction solution;
276:   TSForcingFunction  forcing;

278:   PetscErrorCode (*destroy)(DMTS);
279:   PetscErrorCode (*duplicate)(DMTS,DMTS);
280: };

282: struct _p_DMTS {
283:   PETSCHEADER(struct _DMTSOps);
284:   void *rhsfunctionctx;
285:   void *rhsjacobianctx;

287:   void *ifunctionctx;
288:   void *ijacobianctx;

290:   void *i2functionctx;
291:   void *i2jacobianctx;

293:   void *solutionctx;
294:   void *forcingctx;

296:   void *data;

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

308: PETSC_EXTERN PetscErrorCode DMGetDMTS(DM,DMTS*);
309: PETSC_EXTERN PetscErrorCode DMGetDMTSWrite(DM,DMTS*);
310: PETSC_EXTERN PetscErrorCode DMCopyDMTS(DM,DM);
311: PETSC_EXTERN PetscErrorCode DMTSView(DMTS,PetscViewer);
312: PETSC_EXTERN PetscErrorCode DMTSLoad(DMTS,PetscViewer);
313: PETSC_EXTERN PetscErrorCode DMTSCopy(DMTS,DMTS);

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

317: struct _n_TSEvent {
318:   PetscScalar    *fvalue;          /* value of event function at the end of the step*/
319:   PetscScalar    *fvalue_prev;     /* value of event function at start of the step (left end-point of event interval) */
320:   PetscReal       ptime_prev;      /* time at step start (left end-point of event interval) */
321:   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) */
322:   PetscReal       ptime_right;     /* time on the right end-point of the event interval */
323:   PetscScalar    *fvalue_right;    /* value of event function at the right end-point of the event interval */
324:   PetscInt       *side;            /* Used for detecting repetition of end-point, -1 => left, +1 => right */
325:   PetscReal       timestep_prev;   /* previous time step */
326:   PetscReal       timestep_orig;   /* initial time step */
327:   PetscBool      *zerocrossing;    /* Flag to signal zero crossing detection */
328:   PetscErrorCode  (*eventhandler)(TS,PetscReal,Vec,PetscScalar*,void*); /* User event handler function */
329:   PetscErrorCode  (*postevent)(TS,PetscInt,PetscInt[],PetscReal,Vec,PetscBool,void*); /* User post event function */
330:   void           *ctx;              /* User context for event handler and post even functions */
331:   PetscInt       *direction;        /* Zero crossing direction: 1 -> Going positive, -1 -> Going negative, 0 -> Any */
332:   PetscBool      *terminate;        /* 1 -> Terminate time stepping, 0 -> continue */
333:   PetscInt        nevents;          /* Number of events to handle */
334:   PetscInt        nevents_zero;     /* Number of event zero detected */
335:   PetscInt       *events_zero;      /* List of events that have reached zero */
336:   PetscReal      *vtol;             /* Vector tolerances for event zero check */
337:   TSEventStatus   status;           /* Event status */
338:   PetscInt        iterctr;          /* Iteration counter */
339:   PetscViewer     monitor;
340:   /* Struct to record the events */
341:   struct {
342:     PetscInt  ctr;        /* recorder counter */
343:     PetscReal *time;      /* Event times */
344:     PetscInt  *stepnum;   /* Step numbers */
345:     PetscInt  *nevents;   /* Number of events occuring at the event times */
346:     PetscInt  **eventidx; /* Local indices of the events in the event list */
347:   } recorder;
348:   PetscInt  recsize; /* Size of recorder stack */
349:   PetscInt  refct; /* reference count */
350: };

352: PETSC_EXTERN PetscErrorCode TSEventInitialize(TSEvent,TS,PetscReal,Vec);
353: PETSC_EXTERN PetscErrorCode TSEventDestroy(TSEvent*);
354: PETSC_EXTERN PetscErrorCode TSEventHandler(TS);
355: PETSC_EXTERN PetscErrorCode TSAdjointEventHandler(TS);

357: PETSC_EXTERN PetscLogEvent TS_AdjointStep;
358: PETSC_EXTERN PetscLogEvent TS_Step;
359: PETSC_EXTERN PetscLogEvent TS_PseudoComputeTimeStep;
360: PETSC_EXTERN PetscLogEvent TS_FunctionEval;
361: PETSC_EXTERN PetscLogEvent TS_JacobianEval;
362: PETSC_EXTERN PetscLogEvent TS_ForwardStep;

364: typedef enum {TS_STEP_INCOMPLETE, /* vec_sol, ptime, etc point to beginning of step */
365:               TS_STEP_PENDING,    /* vec_sol advanced, but step has not been accepted yet */
366:               TS_STEP_COMPLETE    /* step accepted and ptime, steps, etc have been advanced */
367: } TSStepStatus;

369: struct _n_TSMonitorLGCtx {
370:   PetscDrawLG    lg;
371:   PetscBool      semilogy;
372:   PetscInt       howoften;  /* when > 0 uses step % howoften, when negative only final solution plotted */
373:   PetscInt       ksp_its,snes_its;
374:   char           **names;
375:   char           **displaynames;
376:   PetscInt       ndisplayvariables;
377:   PetscInt       *displayvariables;
378:   PetscReal      *displayvalues;
379:   PetscErrorCode (*transform)(void*,Vec,Vec*);
380:   PetscErrorCode (*transformdestroy)(void*);
381:   void           *transformctx;
382: };

384: struct _n_TSMonitorEnvelopeCtx {
385:   Vec max,min;
386: };

388: /*
389:     Checks if the user provide a TSSetIFunction() but an explicit method is called; generate an error in that case
390: */
391: PETSC_STATIC_INLINE PetscErrorCode TSCheckImplicitTerm(TS ts)
392: {
393:   TSIFunction      ifunction;
394:   DM               dm;
395:   PetscErrorCode   ierr;

398:   TSGetDM(ts,&dm);
399:   DMTSGetIFunction(dm,&ifunction,NULL);
400:   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()");
401:   return(0);
402: }

404: PETSC_EXTERN PetscLogEvent TSTrajectory_Set;
405: PETSC_EXTERN PetscLogEvent TSTrajectory_Get;
406: PETSC_EXTERN PetscLogEvent TSTrajectory_DiskWrite;
407: PETSC_EXTERN PetscLogEvent TSTrajectory_DiskRead;

409: struct _n_TSMonitorDrawCtx {
410:   PetscViewer   viewer;
411:   Vec           initialsolution;
412:   PetscBool     showinitial;
413:   PetscInt      howoften;  /* when > 0 uses step % howoften, when negative only final solution plotted */
414:   PetscBool     showtimestepandtime;
415: };
416: #endif