xref: /petsc/src/ts/trajectory/impls/memory/trajmemory.c (revision f97672e55eacc8688507b9471cd7ec2664d7f203)
1 #include <petsc/private/tsimpl.h>        /*I "petscts.h"  I*/
2 #include <petscsys.h>
3 #if defined(PETSC_HAVE_REVOLVE)
4 #include <revolve_c.h>
5 
6 /* Limit Revolve to 32-bits */
7 #define PETSC_REVOLVE_INT_MAX  2147483647
8 
9 typedef int PetscRevolveInt;
10 
11 static inline PetscErrorCode PetscRevolveIntCast(PetscInt a,PetscRevolveInt *b)
12 {
13   PetscFunctionBegin;
14 #if defined(PETSC_USE_64BIT_INDICES)
15   *b = 0;
16   PetscCheck(a <= PETSC_REVOLVE_INT_MAX,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Parameter is too large for Revolve, which is restricted to 32 bit integers");
17 #endif
18   *b = (PetscRevolveInt)(a);
19   PetscFunctionReturn(0);
20 }
21 #endif
22 #if defined(PETSC_HAVE_CAMS)
23 #include <offline_schedule.h>
24 #endif
25 
26 PetscLogEvent TSTrajectory_DiskWrite, TSTrajectory_DiskRead;
27 static PetscErrorCode TSTrajectorySet_Memory(TSTrajectory,TS,PetscInt,PetscReal,Vec);
28 
29 typedef enum {NONE,TWO_LEVEL_NOREVOLVE,TWO_LEVEL_REVOLVE,TWO_LEVEL_TWO_REVOLVE,REVOLVE_OFFLINE,REVOLVE_ONLINE,REVOLVE_MULTISTAGE,CAMS_OFFLINE} SchedulerType;
30 
31 typedef enum {UNSET=-1,SOLUTIONONLY=0,STAGESONLY=1,SOLUTION_STAGES=2} CheckpointType;
32 
33 typedef enum {TJ_REVOLVE, TJ_CAMS, TJ_PETSC} TSTrajectoryMemoryType;
34 static const char *const TSTrajectoryMemoryTypes[] = {"REVOLVE","CAMS","PETSC","TSTrajectoryMemoryType","TJ_",NULL};
35 
36 #define HaveSolution(m) ((m) == SOLUTIONONLY || (m) == SOLUTION_STAGES)
37 #define HaveStages(m)   ((m) == STAGESONLY || (m) == SOLUTION_STAGES)
38 
39 typedef struct _StackElement {
40   PetscInt       stepnum;
41   Vec            X;
42   Vec            *Y;
43   PetscReal      time;
44   PetscReal      timeprev; /* for no solution_only mode */
45   PetscReal      timenext; /* for solution_only mode */
46   CheckpointType cptype;
47 } *StackElement;
48 
49 #if defined(PETSC_HAVE_REVOLVE)
50 typedef struct _RevolveCTX {
51   PetscBool       reverseonestep;
52   PetscRevolveInt where;
53   PetscRevolveInt snaps_in;
54   PetscRevolveInt stepsleft;
55   PetscRevolveInt check;
56   PetscRevolveInt oldcapo;
57   PetscRevolveInt capo;
58   PetscRevolveInt fine;
59   PetscRevolveInt info;
60 } RevolveCTX;
61 #endif
62 
63 #if defined(PETSC_HAVE_CAMS)
64 typedef struct _CAMSCTX {
65   PetscInt lastcheckpointstep;
66   PetscInt lastcheckpointtype;
67   PetscInt num_units_avail;
68   PetscInt endstep;
69   PetscInt num_stages;
70   PetscInt nextcheckpointstep;
71   PetscInt nextcheckpointtype; /* (0) solution only (1) stages (2) solution+stages */
72   PetscInt info;
73 } CAMSCTX;
74 #endif
75 
76 typedef struct _Stack {
77   PetscInt      stacksize;
78   PetscInt      top;
79   StackElement  *container;
80   PetscInt      nallocated;
81   PetscInt      numY;
82   PetscBool     solution_only;
83   PetscBool     use_dram;
84 } Stack;
85 
86 typedef struct _DiskStack {
87   PetscInt  stacksize;
88   PetscInt  top;
89   PetscInt  *container;
90 } DiskStack;
91 
92 typedef struct _TJScheduler {
93   SchedulerType stype;
94   TSTrajectoryMemoryType tj_memory_type;
95 #if defined(PETSC_HAVE_REVOLVE)
96   RevolveCTX    *rctx,*rctx2;
97   PetscBool     use_online;
98   PetscInt      store_stride;
99 #endif
100 #if defined(PETSC_HAVE_CAMS)
101   CAMSCTX       *actx;
102 #endif
103   PetscBool     recompute;
104   PetscBool     skip_trajectory;
105   PetscBool     save_stack;
106   PetscInt      max_units_ram;  /* maximum checkpointing units in RAM */
107   PetscInt      max_units_disk; /* maximum checkpointing units on disk */
108   PetscInt      max_cps_ram;    /* maximum checkpoints in RAM */
109   PetscInt      max_cps_disk;   /* maximum checkpoints on disk */
110   PetscInt      stride;
111   PetscInt      total_steps;    /* total number of steps */
112   Stack         stack;
113   DiskStack     diskstack;
114   PetscViewer   viewer;
115 } TJScheduler;
116 
117 static PetscErrorCode TurnForwardWithStepsize(TS ts,PetscReal nextstepsize)
118 {
119   PetscFunctionBegin;
120   /* reverse the direction */
121   PetscCall(TSSetTimeStep(ts,nextstepsize));
122   PetscFunctionReturn(0);
123 }
124 
125 static PetscErrorCode TurnForward(TS ts)
126 {
127   PetscReal      stepsize;
128 
129   PetscFunctionBegin;
130   /* reverse the direction */
131   PetscCall(TSGetTimeStep(ts,&stepsize));
132   PetscCall(TSSetTimeStep(ts,-stepsize));
133   PetscFunctionReturn(0);
134 }
135 
136 static PetscErrorCode TurnBackward(TS ts)
137 {
138   PetscReal      stepsize;
139 
140   PetscFunctionBegin;
141   if (!ts->trajectory->adjoint_solve_mode) PetscFunctionReturn(0);
142   /* reverse the direction */
143   stepsize = ts->ptime_prev-ts->ptime;
144   PetscCall(TSSetTimeStep(ts,stepsize));
145   PetscFunctionReturn(0);
146 }
147 
148 static PetscErrorCode ElementCreate(TS ts,CheckpointType cptype,Stack *stack,StackElement *e)
149 {
150   Vec            X;
151   Vec            *Y;
152 
153   PetscFunctionBegin;
154   if (stack->top < stack->stacksize-1 && stack->container[stack->top+1]) {
155     *e = stack->container[stack->top+1];
156     if (HaveSolution(cptype) && !(*e)->X) {
157       PetscCall(TSGetSolution(ts,&X));
158       PetscCall(VecDuplicate(X,&(*e)->X));
159     }
160     if (cptype==1 && (*e)->X) {
161       PetscCall(VecDestroy(&(*e)->X));
162     }
163     if (HaveStages(cptype) && !(*e)->Y) {
164       PetscCall(TSGetStages(ts,&stack->numY,&Y));
165       if (stack->numY) {
166         PetscCall(VecDuplicateVecs(Y[0],stack->numY,&(*e)->Y));
167       }
168     }
169     if (cptype==0 && (*e)->Y) {
170       PetscCall(VecDestroyVecs(stack->numY,&(*e)->Y));
171     }
172     (*e)->cptype = cptype;
173     PetscFunctionReturn(0);
174   }
175   if (stack->use_dram) {
176     PetscCall(PetscMallocSetDRAM());
177   }
178   PetscCall(PetscNew(e));
179   if (HaveSolution(cptype)) {
180     PetscCall(TSGetSolution(ts,&X));
181     PetscCall(VecDuplicate(X,&(*e)->X));
182   }
183   if (HaveStages(cptype)) {
184     PetscCall(TSGetStages(ts,&stack->numY,&Y));
185     if (stack->numY) {
186       PetscCall(VecDuplicateVecs(Y[0],stack->numY,&(*e)->Y));
187     }
188   }
189   if (stack->use_dram) {
190     PetscCall(PetscMallocResetDRAM());
191   }
192   stack->nallocated++;
193   (*e)->cptype = cptype;
194   PetscFunctionReturn(0);
195 }
196 
197 static PetscErrorCode ElementSet(TS ts, Stack *stack, StackElement *e, PetscInt stepnum, PetscReal time, Vec X)
198 {
199   Vec            *Y;
200   PetscInt       i;
201   PetscReal      timeprev;
202 
203   PetscFunctionBegin;
204   if (HaveSolution((*e)->cptype)) {
205     PetscCall(VecCopy(X,(*e)->X));
206   }
207   if (HaveStages((*e)->cptype)) {
208     PetscCall(TSGetStages(ts,&stack->numY,&Y));
209     for (i=0;i<stack->numY;i++) {
210       PetscCall(VecCopy(Y[i],(*e)->Y[i]));
211     }
212   }
213   (*e)->stepnum = stepnum;
214   (*e)->time    = time;
215   /* for consistency */
216   if (stepnum == 0) {
217     (*e)->timeprev = (*e)->time - ts->time_step;
218   } else {
219     PetscCall(TSGetPrevTime(ts,&timeprev));
220     (*e)->timeprev = timeprev;
221   }
222   PetscFunctionReturn(0);
223 }
224 
225 static PetscErrorCode ElementDestroy(Stack *stack,StackElement e)
226 {
227   PetscFunctionBegin;
228   if (stack->use_dram) {
229     PetscCall(PetscMallocSetDRAM());
230   }
231   PetscCall(VecDestroy(&e->X));
232   if (e->Y) {
233     PetscCall(VecDestroyVecs(stack->numY,&e->Y));
234   }
235   PetscCall(PetscFree(e));
236   if (stack->use_dram) {
237     PetscCall(PetscMallocResetDRAM());
238   }
239   stack->nallocated--;
240   PetscFunctionReturn(0);
241 }
242 
243 static PetscErrorCode StackResize(Stack *stack,PetscInt newsize)
244 {
245   StackElement   *newcontainer;
246   PetscInt       i;
247 
248   PetscFunctionBegin;
249   PetscCall(PetscCalloc1(newsize*sizeof(StackElement),&newcontainer));
250   for (i=0;i<stack->stacksize;i++) {
251     newcontainer[i] = stack->container[i];
252   }
253   PetscCall(PetscFree(stack->container));
254   stack->container = newcontainer;
255   stack->stacksize = newsize;
256   PetscFunctionReturn(0);
257 }
258 
259 static PetscErrorCode StackPush(Stack *stack,StackElement e)
260 {
261   PetscFunctionBegin;
262   PetscCheck(stack->top+1 < stack->stacksize,PETSC_COMM_SELF,PETSC_ERR_MEMC,"Maximum stack size (%" PetscInt_FMT ") exceeded",stack->stacksize);
263   stack->container[++stack->top] = e;
264   PetscFunctionReturn(0);
265 }
266 
267 static PetscErrorCode StackPop(Stack *stack,StackElement *e)
268 {
269   PetscFunctionBegin;
270   *e = NULL;
271   PetscCheck(stack->top != -1,PETSC_COMM_SELF,PETSC_ERR_MEMC,"Empty stack");
272   *e = stack->container[stack->top--];
273   PetscFunctionReturn(0);
274 }
275 
276 static PetscErrorCode StackTop(Stack *stack,StackElement *e)
277 {
278   PetscFunctionBegin;
279   *e = stack->container[stack->top];
280   PetscFunctionReturn(0);
281 }
282 
283 static PetscErrorCode StackInit(Stack *stack,PetscInt size,PetscInt ny)
284 {
285   PetscFunctionBegin;
286   stack->top  = -1;
287   stack->numY = ny;
288 
289   if (!stack->container) {
290     PetscCall(PetscCalloc1(size,&stack->container));
291   }
292   PetscFunctionReturn(0);
293 }
294 
295 static PetscErrorCode StackDestroy(Stack *stack)
296 {
297   const PetscInt n = stack->nallocated;
298 
299   PetscFunctionBegin;
300   if (!stack->container) PetscFunctionReturn(0);
301   PetscCheck(stack->top+1 <= n,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Stack size does not match element counter %" PetscInt_FMT,n);
302   for (PetscInt i=0; i<n; i++) PetscCall(ElementDestroy(stack,stack->container[i]));
303   PetscCall(PetscFree(stack->container));
304   PetscFunctionReturn(0);
305 }
306 
307 static PetscErrorCode StackFind(Stack *stack,StackElement *e,PetscInt index)
308 {
309   PetscFunctionBegin;
310   *e = NULL;
311   PetscCheck(index >= 0 && index <= stack->top,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Invalid index %" PetscInt_FMT,index);
312   *e = stack->container[index];
313   PetscFunctionReturn(0);
314 }
315 
316 static PetscErrorCode WriteToDisk(PetscBool stifflyaccurate,PetscInt stepnum,PetscReal time,PetscReal timeprev,Vec X,Vec *Y,PetscInt numY,CheckpointType cptype,PetscViewer viewer)
317 {
318   PetscFunctionBegin;
319   PetscCall(PetscViewerBinaryWrite(viewer,&stepnum,1,PETSC_INT));
320   if (HaveSolution(cptype)) PetscCall(VecView(X,viewer));
321   if (HaveStages(cptype)) {
322     for (PetscInt i=0; i<numY; i++) {
323       /* For stiffly accurate TS methods, the last stage Y[ns-1] is the same as the solution X, thus does not need to be saved again. */
324       if (stifflyaccurate && i == numY-1 && HaveSolution(cptype)) continue;
325       PetscCall(VecView(Y[i],viewer));
326     }
327   }
328   PetscCall(PetscViewerBinaryWrite(viewer,&time,1,PETSC_REAL));
329   PetscCall(PetscViewerBinaryWrite(viewer,&timeprev,1,PETSC_REAL));
330   PetscFunctionReturn(0);
331 }
332 
333 static PetscErrorCode ReadFromDisk(PetscBool stifflyaccurate,PetscInt *stepnum,PetscReal *time,PetscReal *timeprev,Vec X,Vec *Y,PetscInt numY,CheckpointType cptype,PetscViewer viewer)
334 {
335   PetscFunctionBegin;
336   PetscCall(PetscViewerBinaryRead(viewer,stepnum,1,NULL,PETSC_INT));
337   if (HaveSolution(cptype)) PetscCall(VecLoad(X,viewer));
338   if (HaveStages(cptype)) {
339     for (PetscInt i=0; i<numY; i++) {
340       /* For stiffly accurate TS methods, the last stage Y[ns-1] is the same as the solution X, thus does not need to be loaded again. */
341       if (stifflyaccurate && i == numY-1 && HaveSolution(cptype)) continue;
342       PetscCall(VecLoad(Y[i],viewer));
343     }
344   }
345   PetscCall(PetscViewerBinaryRead(viewer,time,1,NULL,PETSC_REAL));
346   PetscCall(PetscViewerBinaryRead(viewer,timeprev,1,NULL,PETSC_REAL));
347   PetscFunctionReturn(0);
348 }
349 
350 static PetscErrorCode StackDumpAll(TSTrajectory tj,TS ts,Stack *stack,PetscInt id)
351 {
352   Vec          *Y;
353   PetscInt      ndumped,cptype_int;
354   StackElement  e     = NULL;
355   TJScheduler  *tjsch = (TJScheduler*)tj->data;
356   char          filename[PETSC_MAX_PATH_LEN];
357   MPI_Comm      comm;
358 
359   PetscFunctionBegin;
360   PetscCall(PetscObjectGetComm((PetscObject)ts,&comm));
361   if (tj->monitor) {
362     PetscCall(PetscViewerASCIIPushTab(tj->monitor));
363     PetscCall(PetscViewerASCIIPrintf(tj->monitor,"Dump stack id %" PetscInt_FMT " to file\n",id));
364     PetscCall(PetscViewerASCIIPopTab(tj->monitor));
365   }
366   PetscCall(PetscSNPrintf(filename,sizeof(filename),"%s/TS-STACK%06" PetscInt_FMT ".bin",tj->dirname,id));
367   PetscCall(PetscViewerFileSetName(tjsch->viewer,filename));
368   PetscCall(PetscViewerSetUp(tjsch->viewer));
369   ndumped = stack->top+1;
370   PetscCall(PetscViewerBinaryWrite(tjsch->viewer,&ndumped,1,PETSC_INT));
371   for (PetscInt i=0;i<ndumped;i++) {
372     e = stack->container[i];
373     cptype_int = (PetscInt)e->cptype;
374     PetscCall(PetscViewerBinaryWrite(tjsch->viewer,&cptype_int,1,PETSC_INT));
375     PetscCall(PetscLogEventBegin(TSTrajectory_DiskWrite,tj,ts,0,0));
376     PetscCall(WriteToDisk(ts->stifflyaccurate,e->stepnum,e->time,e->timeprev,e->X,e->Y,stack->numY,e->cptype,tjsch->viewer));
377     PetscCall(PetscLogEventEnd(TSTrajectory_DiskWrite,tj,ts,0,0));
378     ts->trajectory->diskwrites++;
379     PetscCall(StackPop(stack,&e));
380   }
381   /* save the last step for restart, the last step is in memory when using single level schemes, but not necessarily the case for multi level schemes */
382   PetscCall(TSGetStages(ts,&stack->numY,&Y));
383   PetscCall(PetscLogEventBegin(TSTrajectory_DiskWrite,tj,ts,0,0));
384   PetscCall(WriteToDisk(ts->stifflyaccurate,ts->steps,ts->ptime,ts->ptime_prev,ts->vec_sol,Y,stack->numY,SOLUTION_STAGES,tjsch->viewer));
385   PetscCall(PetscLogEventEnd(TSTrajectory_DiskWrite,tj,ts,0,0));
386   ts->trajectory->diskwrites++;
387   PetscFunctionReturn(0);
388 }
389 
390 static PetscErrorCode StackLoadAll(TSTrajectory tj,TS ts,Stack *stack,PetscInt id)
391 {
392   Vec          *Y;
393   PetscInt      i,nloaded,cptype_int;
394   StackElement  e;
395   PetscViewer   viewer;
396   char          filename[PETSC_MAX_PATH_LEN];
397 
398   PetscFunctionBegin;
399   if (tj->monitor) {
400     PetscCall(PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel));
401     PetscCall(PetscViewerASCIIPrintf(tj->monitor,"Load stack from file\n"));
402     PetscCall(PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel));
403   }
404   PetscCall(PetscSNPrintf(filename,sizeof filename,"%s/TS-STACK%06" PetscInt_FMT ".bin",tj->dirname,id));
405   PetscCall(PetscViewerBinaryOpen(PetscObjectComm((PetscObject)tj),filename,FILE_MODE_READ,&viewer));
406   PetscCall(PetscViewerBinarySetSkipInfo(viewer,PETSC_TRUE));
407   PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_NATIVE));
408   PetscCall(PetscViewerBinaryRead(viewer,&nloaded,1,NULL,PETSC_INT));
409   for (i=0;i<nloaded;i++) {
410     PetscCall(PetscViewerBinaryRead(viewer,&cptype_int,1,NULL,PETSC_INT));
411     PetscCall(ElementCreate(ts,(CheckpointType)cptype_int,stack,&e));
412     PetscCall(StackPush(stack,e));
413     PetscCall(PetscLogEventBegin(TSTrajectory_DiskRead,tj,ts,0,0));
414     PetscCall(ReadFromDisk(ts->stifflyaccurate,&e->stepnum,&e->time,&e->timeprev,e->X,e->Y,stack->numY,e->cptype,viewer));
415     PetscCall(PetscLogEventEnd(TSTrajectory_DiskRead,tj,ts,0,0));
416     ts->trajectory->diskreads++;
417   }
418   /* load the last step into TS */
419   PetscCall(TSGetStages(ts,&stack->numY,&Y));
420   PetscCall(PetscLogEventBegin(TSTrajectory_DiskRead,tj,ts,0,0));
421   PetscCall(ReadFromDisk(ts->stifflyaccurate,&ts->steps,&ts->ptime,&ts->ptime_prev,ts->vec_sol,Y,stack->numY,SOLUTION_STAGES,viewer));
422   PetscCall(PetscLogEventEnd(TSTrajectory_DiskRead,tj,ts,0,0));
423   ts->trajectory->diskreads++;
424   PetscCall(TurnBackward(ts));
425   PetscCall(PetscViewerDestroy(&viewer));
426   PetscFunctionReturn(0);
427 }
428 
429 #if defined(PETSC_HAVE_REVOLVE)
430 static PetscErrorCode StackLoadLast(TSTrajectory tj,TS ts,Stack *stack,PetscInt id)
431 {
432   Vec         *Y;
433   PetscInt     size;
434   PetscViewer  viewer;
435   char         filename[PETSC_MAX_PATH_LEN];
436 #if defined(PETSC_HAVE_MPIIO)
437   PetscBool    usempiio;
438 #endif
439   int          fd;
440   off_t        off,offset;
441 
442   PetscFunctionBegin;
443   if (tj->monitor) {
444     PetscCall(PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel));
445     PetscCall(PetscViewerASCIIPrintf(tj->monitor,"Load last stack element from file\n"));
446     PetscCall(PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel));
447   }
448   PetscCall(TSGetStages(ts,&stack->numY,&Y));
449   PetscCall(VecGetSize(Y[0],&size));
450   /* VecView writes to file two extra int's for class id and number of rows */
451   off  = -((stack->solution_only?0:stack->numY)+1)*(size*PETSC_BINARY_SCALAR_SIZE+2*PETSC_BINARY_INT_SIZE)-PETSC_BINARY_INT_SIZE-2*PETSC_BINARY_SCALAR_SIZE;
452 
453   PetscCall(PetscSNPrintf(filename,sizeof filename,"%s/TS-STACK%06" PetscInt_FMT ".bin",tj->dirname,id));
454   PetscCall(PetscViewerBinaryOpen(PetscObjectComm((PetscObject)tj),filename,FILE_MODE_READ,&viewer));
455   PetscCall(PetscViewerBinarySetSkipInfo(viewer,PETSC_TRUE));
456   PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_NATIVE));
457 #if defined(PETSC_HAVE_MPIIO)
458   PetscCall(PetscViewerBinaryGetUseMPIIO(viewer,&usempiio));
459   if (usempiio) {
460     PetscCall(PetscViewerBinaryGetMPIIODescriptor(viewer,(MPI_File*)&fd));
461     PetscCall(PetscBinarySynchronizedSeek(PetscObjectComm((PetscObject)tj),fd,off,PETSC_BINARY_SEEK_END,&offset));
462   } else {
463 #endif
464     PetscCall(PetscViewerBinaryGetDescriptor(viewer,&fd));
465     PetscCall(PetscBinarySeek(fd,off,PETSC_BINARY_SEEK_END,&offset));
466 #if defined(PETSC_HAVE_MPIIO)
467   }
468 #endif
469   /* load the last step into TS */
470   PetscCall(PetscLogEventBegin(TSTrajectory_DiskRead,tj,ts,0,0));
471   PetscCall(ReadFromDisk(ts->stifflyaccurate,&ts->steps,&ts->ptime,&ts->ptime_prev,ts->vec_sol,Y,stack->numY,SOLUTION_STAGES,viewer));
472   PetscCall(PetscLogEventEnd(TSTrajectory_DiskRead,tj,ts,0,0));
473   ts->trajectory->diskreads++;
474   PetscCall(PetscViewerDestroy(&viewer));
475   PetscCall(TurnBackward(ts));
476   PetscFunctionReturn(0);
477 }
478 #endif
479 
480 static PetscErrorCode DumpSingle(TSTrajectory tj,TS ts,Stack *stack,PetscInt id)
481 {
482   Vec            *Y;
483   PetscInt       stepnum;
484   TJScheduler    *tjsch = (TJScheduler*)tj->data;
485   char           filename[PETSC_MAX_PATH_LEN];
486   MPI_Comm       comm;
487 
488   PetscFunctionBegin;
489   PetscCall(PetscObjectGetComm((PetscObject)ts,&comm));
490   if (tj->monitor) {
491     PetscCall(PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel));
492     PetscCall(PetscViewerASCIIPrintf(tj->monitor,"Dump a single point from file\n"));
493     PetscCall(PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel));
494   }
495   PetscCall(TSGetStepNumber(ts,&stepnum));
496   PetscCall(PetscSNPrintf(filename,sizeof(filename),"%s/TS-CPS%06" PetscInt_FMT ".bin",tj->dirname,id));
497   PetscCall(PetscViewerFileSetName(tjsch->viewer,filename));
498   PetscCall(PetscViewerSetUp(tjsch->viewer));
499 
500   PetscCall(TSGetStages(ts,&stack->numY,&Y));
501   PetscCall(PetscLogEventBegin(TSTrajectory_DiskWrite,tj,ts,0,0));
502   PetscCall(WriteToDisk(ts->stifflyaccurate,stepnum,ts->ptime,ts->ptime_prev,ts->vec_sol,Y,stack->numY,SOLUTION_STAGES,tjsch->viewer));
503   PetscCall(PetscLogEventEnd(TSTrajectory_DiskWrite,tj,ts,0,0));
504   ts->trajectory->diskwrites++;
505   PetscFunctionReturn(0);
506 }
507 
508 static PetscErrorCode LoadSingle(TSTrajectory tj,TS ts,Stack *stack,PetscInt id)
509 {
510   Vec            *Y;
511   PetscViewer    viewer;
512   char           filename[PETSC_MAX_PATH_LEN];
513 
514   PetscFunctionBegin;
515   if (tj->monitor) {
516     PetscCall(PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel));
517     PetscCall(PetscViewerASCIIPrintf(tj->monitor,"Load a single point from file\n"));
518     PetscCall(PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel));
519   }
520   PetscCall(PetscSNPrintf(filename,sizeof filename,"%s/TS-CPS%06" PetscInt_FMT ".bin",tj->dirname,id));
521   PetscCall(PetscViewerBinaryOpen(PetscObjectComm((PetscObject)tj),filename,FILE_MODE_READ,&viewer));
522   PetscCall(PetscViewerBinarySetSkipInfo(viewer,PETSC_TRUE));
523   PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_NATIVE));
524   PetscCall(TSGetStages(ts,&stack->numY,&Y));
525   PetscCall(PetscLogEventBegin(TSTrajectory_DiskRead,tj,ts,0,0));
526   PetscCall(ReadFromDisk(ts->stifflyaccurate,&ts->steps,&ts->ptime,&ts->ptime_prev,ts->vec_sol,Y,stack->numY,SOLUTION_STAGES,viewer));
527   PetscCall(PetscLogEventEnd(TSTrajectory_DiskRead,tj,ts,0,0));
528   ts->trajectory->diskreads++;
529   PetscCall(PetscViewerDestroy(&viewer));
530   PetscFunctionReturn(0);
531 }
532 
533 static PetscErrorCode UpdateTS(TS ts,Stack *stack,StackElement e,PetscInt stepnum,PetscBool adjoint_mode)
534 {
535   Vec            *Y;
536   PetscInt       i;
537 
538   PetscFunctionBegin;
539   /* In adjoint mode we do not need to copy solution if the stepnum is the same */
540   if (!adjoint_mode || (HaveSolution(e->cptype) && e->stepnum!=stepnum)) {
541     PetscCall(VecCopy(e->X,ts->vec_sol));
542   }
543   if (HaveStages(e->cptype)) {
544     PetscCall(TSGetStages(ts,&stack->numY,&Y));
545     if (e->stepnum && e->stepnum==stepnum) {
546       for (i=0;i<stack->numY;i++) {
547         PetscCall(VecCopy(e->Y[i],Y[i]));
548       }
549     } else if (ts->stifflyaccurate) {
550       PetscCall(VecCopy(e->Y[stack->numY-1],ts->vec_sol));
551     }
552   }
553   if (adjoint_mode) {
554     PetscCall(TSSetTimeStep(ts,e->timeprev-e->time)); /* stepsize will be negative */
555   } else {
556     PetscCall(TSSetTimeStep(ts,e->time-e->timeprev)); /* stepsize will be positive */
557   }
558   ts->ptime      = e->time;
559   ts->ptime_prev = e->timeprev;
560   PetscFunctionReturn(0);
561 }
562 
563 static PetscErrorCode ReCompute(TS ts,TJScheduler *tjsch,PetscInt stepnumbegin,PetscInt stepnumend)
564 {
565   Stack          *stack = &tjsch->stack;
566   PetscInt       i;
567 
568   PetscFunctionBegin;
569   tjsch->recompute = PETSC_TRUE; /* hints TSTrajectorySet() that it is in recompute mode */
570   PetscCall(TSSetStepNumber(ts,stepnumbegin));/* global step number */
571   for (i=stepnumbegin;i<stepnumend;i++) { /* assume fixed step size */
572     if (stack->solution_only && !tjsch->skip_trajectory) { /* revolve online need this */
573       /* don't use the public interface as it will update the TSHistory: this need a better fix */
574       PetscCall(TSTrajectorySet_Memory(ts->trajectory,ts,ts->steps,ts->ptime,ts->vec_sol));
575     }
576     PetscCall(TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol));
577     PetscCall(TSStep(ts));
578     if (!stack->solution_only && !tjsch->skip_trajectory) {
579       /* don't use the public interface as it will update the TSHistory: this need a better fix */
580       PetscCall(TSTrajectorySet_Memory(ts->trajectory,ts,ts->steps,ts->ptime,ts->vec_sol));
581     }
582     PetscCall(TSEventHandler(ts));
583     if (!ts->steprollback) {
584       PetscCall(TSPostStep(ts));
585     }
586   }
587   PetscCall(TurnBackward(ts));
588   ts->trajectory->recomps += stepnumend-stepnumbegin; /* recomputation counter */
589   PetscCall(TSSetStepNumber(ts,stepnumend));
590   tjsch->recompute = PETSC_FALSE; /* reset the flag for recompute mode */
591   PetscFunctionReturn(0);
592 }
593 
594 static PetscErrorCode TopLevelStore(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscInt localstepnum,PetscInt laststridesize,PetscBool *done)
595 {
596   Stack          *stack = &tjsch->stack;
597   DiskStack      *diskstack = &tjsch->diskstack;
598   PetscInt       stridenum;
599 
600   PetscFunctionBegin;
601   *done = PETSC_FALSE;
602   stridenum    = stepnum/tjsch->stride;
603   /* make sure saved checkpoint id starts from 1
604      skip last stride when using stridenum+1
605      skip first stride when using stridenum */
606   if (stack->solution_only) {
607     if (tjsch->save_stack) {
608       if (localstepnum == tjsch->stride-1 && stepnum < tjsch->total_steps-laststridesize) { /* current step will be saved without going through stack */
609         PetscCall(StackDumpAll(tj,ts,stack,stridenum+1));
610         if (tjsch->stype == TWO_LEVEL_TWO_REVOLVE) diskstack->container[++diskstack->top] = stridenum+1;
611         *done = PETSC_TRUE;
612       }
613     } else {
614       if (localstepnum == 0 && stepnum < tjsch->total_steps-laststridesize) {
615         PetscCall(DumpSingle(tj,ts,stack,stridenum+1));
616         if (tjsch->stype == TWO_LEVEL_TWO_REVOLVE) diskstack->container[++diskstack->top] = stridenum+1;
617         *done = PETSC_TRUE;
618       }
619     }
620   } else {
621     if (tjsch->save_stack) {
622       if (localstepnum == 0 && stepnum < tjsch->total_steps && stepnum != 0) { /* skip the first stride */
623         PetscCall(StackDumpAll(tj,ts,stack,stridenum));
624         if (tjsch->stype == TWO_LEVEL_TWO_REVOLVE) diskstack->container[++diskstack->top] = stridenum;
625         *done = PETSC_TRUE;
626       }
627     } else {
628       if (localstepnum == 1 && stepnum < tjsch->total_steps-laststridesize) {
629         PetscCall(DumpSingle(tj,ts,stack,stridenum+1));
630         if (tjsch->stype == TWO_LEVEL_TWO_REVOLVE) diskstack->container[++diskstack->top] = stridenum+1;
631         *done = PETSC_TRUE;
632       }
633     }
634   }
635   PetscFunctionReturn(0);
636 }
637 
638 static PetscErrorCode TSTrajectoryMemorySet_N(TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X)
639 {
640   Stack          *stack = &tjsch->stack;
641   StackElement   e;
642   CheckpointType cptype;
643 
644   PetscFunctionBegin;
645   /* skip the last step */
646   if (ts->reason) { /* only affect the forward run */
647     /* update total_steps in the end of forward run */
648     if (stepnum != tjsch->total_steps) tjsch->total_steps = stepnum;
649     if (stack->solution_only) {
650       /* get rid of the solution at second last step */
651       PetscCall(StackPop(stack,&e));
652     }
653     PetscFunctionReturn(0);
654   }
655   /*  do not save trajectory at the recompute stage for solution_only mode */
656   if (stack->solution_only && tjsch->recompute) PetscFunctionReturn(0);
657   /* skip the first step for no_solution_only mode */
658   if (!stack->solution_only && stepnum == 0) PetscFunctionReturn(0);
659 
660   /* resize the stack */
661   if (stack->top+1 == stack->stacksize) {
662     PetscCall(StackResize(stack,2*stack->stacksize));
663   }
664   /* update timenext for the previous step; necessary for step adaptivity */
665   if (stack->top > -1) {
666     PetscCall(StackTop(stack,&e));
667     e->timenext = ts->ptime;
668   }
669   if (stepnum < stack->top) {
670     SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_MEMC,"Illegal modification of a non-top stack element");
671   }
672   cptype = stack->solution_only ? SOLUTIONONLY : STAGESONLY;
673   PetscCall(ElementCreate(ts,cptype,stack,&e));
674   PetscCall(ElementSet(ts,stack,&e,stepnum,time,X));
675   PetscCall(StackPush(stack,e));
676   PetscFunctionReturn(0);
677 }
678 
679 static PetscErrorCode TSTrajectoryMemorySet_N_2(TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X)
680 {
681   Stack          *stack = &tjsch->stack;
682   StackElement   e;
683   CheckpointType cptype;
684 
685   PetscFunctionBegin;
686   if (stack->top+1 == stack->stacksize) {
687     PetscCall(StackResize(stack,2*stack->stacksize));
688   }
689   /* update timenext for the previous step; necessary for step adaptivity */
690   if (stack->top > -1) {
691     PetscCall(StackTop(stack,&e));
692     e->timenext = ts->ptime;
693   }
694   PetscCheck(stepnum >= stack->top,PetscObjectComm((PetscObject)ts),PETSC_ERR_MEMC,"Illegal modification of a non-top stack element");
695   cptype = stack->solution_only ? SOLUTIONONLY : SOLUTION_STAGES; /* Always include solution in a checkpoint in non-adjoint mode */
696   PetscCall(ElementCreate(ts,cptype,stack,&e));
697   PetscCall(ElementSet(ts,stack,&e,stepnum,time,X));
698   PetscCall(StackPush(stack,e));
699   PetscFunctionReturn(0);
700 }
701 
702 static PetscErrorCode TSTrajectoryMemoryGet_N(TS ts,TJScheduler *tjsch,PetscInt stepnum)
703 {
704   Stack          *stack = &tjsch->stack;
705   StackElement   e;
706   PetscInt       ns;
707 
708   PetscFunctionBegin;
709   /* If TSTrajectoryGet() is called after TSAdjointSolve() converges (e.g. outside the while loop in TSAdjointSolve()), skip getting the checkpoint. */
710   if (ts->reason) PetscFunctionReturn(0);
711   if (stepnum == tjsch->total_steps) {
712     PetscCall(TurnBackward(ts));
713     PetscFunctionReturn(0);
714   }
715   /* restore a checkpoint */
716   PetscCall(StackTop(stack,&e));
717   PetscCall(UpdateTS(ts,stack,e,stepnum,PETSC_TRUE));
718   PetscCall(TSGetStages(ts,&ns,PETSC_IGNORE));
719   if (stack->solution_only && ns) { /* recompute one step */
720     PetscCall(TurnForwardWithStepsize(ts,e->timenext-e->time));
721     PetscCall(ReCompute(ts,tjsch,e->stepnum,stepnum));
722   }
723   PetscCall(StackPop(stack,&e));
724   PetscFunctionReturn(0);
725 }
726 
727 static PetscErrorCode TSTrajectoryMemoryGet_N_2(TS ts,TJScheduler *tjsch,PetscInt stepnum)
728 {
729   Stack        *stack = &tjsch->stack;
730   StackElement  e     = NULL;
731 
732   PetscFunctionBegin;
733   PetscCall(StackFind(stack,&e,stepnum));
734   PetscCheck(stepnum == e->stepnum,PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Inconsistent steps! %" PetscInt_FMT " != %" PetscInt_FMT,stepnum,e->stepnum);
735   PetscCall(UpdateTS(ts,stack,e,stepnum,PETSC_FALSE));
736   PetscFunctionReturn(0);
737 }
738 
739 static PetscErrorCode TSTrajectoryMemorySet_TLNR(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X)
740 {
741   Stack          *stack = &tjsch->stack;
742   PetscInt       localstepnum,laststridesize;
743   StackElement   e;
744   PetscBool      done;
745   CheckpointType cptype;
746 
747   PetscFunctionBegin;
748   if (!stack->solution_only && stepnum == 0) PetscFunctionReturn(0);
749   if (stack->solution_only && stepnum == tjsch->total_steps) PetscFunctionReturn(0);
750   if (tjsch->save_stack && tjsch->recompute) PetscFunctionReturn(0);
751 
752   localstepnum = stepnum%tjsch->stride;
753   /* (stridesize-1) checkpoints are saved in each stride; an extra point is added by StackDumpAll() */
754   laststridesize = tjsch->total_steps%tjsch->stride;
755   if (!laststridesize) laststridesize = tjsch->stride;
756 
757   if (!tjsch->recompute) {
758     PetscCall(TopLevelStore(tj,ts,tjsch,stepnum,localstepnum,laststridesize,&done));
759     if (!tjsch->save_stack && stepnum < tjsch->total_steps-laststridesize) PetscFunctionReturn(0);
760   }
761   if (!stack->solution_only && localstepnum == 0) PetscFunctionReturn(0); /* skip last point in each stride at recompute stage or last stride */
762   if (stack->solution_only && localstepnum == tjsch->stride-1) PetscFunctionReturn(0); /* skip last step in each stride at recompute stage or last stride */
763 
764   cptype = stack->solution_only ? SOLUTIONONLY : STAGESONLY;
765   PetscCall(ElementCreate(ts,cptype,stack,&e));
766   PetscCall(ElementSet(ts,stack,&e,stepnum,time,X));
767   PetscCall(StackPush(stack,e));
768   PetscFunctionReturn(0);
769 }
770 
771 static PetscErrorCode TSTrajectoryMemoryGet_TLNR(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum)
772 {
773   Stack          *stack = &tjsch->stack;
774   PetscInt       id,localstepnum,laststridesize;
775   StackElement   e;
776 
777   PetscFunctionBegin;
778   if (stepnum == tjsch->total_steps) {
779     PetscCall(TurnBackward(ts));
780     PetscFunctionReturn(0);
781   }
782 
783   localstepnum = stepnum%tjsch->stride;
784   laststridesize = tjsch->total_steps%tjsch->stride;
785   if (!laststridesize) laststridesize = tjsch->stride;
786   if (stack->solution_only) {
787     /* fill stack with info */
788     if (localstepnum == 0 && tjsch->total_steps-stepnum >= laststridesize) {
789       id = stepnum/tjsch->stride;
790       if (tjsch->save_stack) {
791         PetscCall(StackLoadAll(tj,ts,stack,id));
792         tjsch->skip_trajectory = PETSC_TRUE;
793         PetscCall(TurnForward(ts));
794         PetscCall(ReCompute(ts,tjsch,id*tjsch->stride-1,id*tjsch->stride));
795         tjsch->skip_trajectory = PETSC_FALSE;
796       } else {
797         PetscCall(LoadSingle(tj,ts,stack,id));
798         PetscCall(TurnForward(ts));
799         PetscCall(ReCompute(ts,tjsch,(id-1)*tjsch->stride,id*tjsch->stride));
800       }
801       PetscFunctionReturn(0);
802     }
803     /* restore a checkpoint */
804     PetscCall(StackPop(stack,&e));
805     PetscCall(UpdateTS(ts,stack,e,stepnum,PETSC_TRUE));
806     tjsch->skip_trajectory = PETSC_TRUE;
807     PetscCall(TurnForward(ts));
808     PetscCall(ReCompute(ts,tjsch,e->stepnum,stepnum));
809     tjsch->skip_trajectory = PETSC_FALSE;
810   } else {
811     CheckpointType cptype = STAGESONLY;
812     /* fill stack with info */
813     if (localstepnum == 0 && tjsch->total_steps-stepnum >= laststridesize) {
814       id = stepnum/tjsch->stride;
815       if (tjsch->save_stack) {
816         PetscCall(StackLoadAll(tj,ts,stack,id));
817       } else {
818         PetscCall(LoadSingle(tj,ts,stack,id));
819         PetscCall(ElementCreate(ts,cptype,stack,&e));
820         PetscCall(ElementSet(ts,stack,&e,(id-1)*tjsch->stride+1,ts->ptime,ts->vec_sol));
821         PetscCall(StackPush(stack,e));
822         PetscCall(TurnForward(ts));
823         PetscCall(ReCompute(ts,tjsch,e->stepnum,id*tjsch->stride));
824       }
825       PetscFunctionReturn(0);
826     }
827     /* restore a checkpoint */
828     PetscCall(StackPop(stack,&e));
829     PetscCall(UpdateTS(ts,stack,e,stepnum,PETSC_TRUE));
830   }
831   PetscFunctionReturn(0);
832 }
833 
834 #if defined(PETSC_HAVE_REVOLVE)
835 static PetscErrorCode printwhattodo(PetscViewer viewer,PetscRevolveInt whattodo,RevolveCTX *rctx,PetscRevolveInt shift)
836 {
837   PetscFunctionBegin;
838   if (!viewer) PetscFunctionReturn(0);
839 
840   switch(whattodo) {
841     case 1:
842       PetscCall(PetscViewerASCIIPrintf(viewer,"Advance from %d to %d\n",rctx->oldcapo+shift,rctx->capo+shift));
843       break;
844     case 2:
845       PetscCall(PetscViewerASCIIPrintf(viewer,"Store in checkpoint number %d (located in RAM)\n",rctx->check));
846       break;
847     case 3:
848       PetscCall(PetscViewerASCIIPrintf(viewer,"First turn: Initialize adjoints and reverse first step\n"));
849       break;
850     case 4:
851       PetscCall(PetscViewerASCIIPrintf(viewer,"Forward and reverse one step\n"));
852       break;
853     case 5:
854       PetscCall(PetscViewerASCIIPrintf(viewer,"Restore in checkpoint number %d (located in RAM)\n",rctx->check));
855       break;
856     case 7:
857       PetscCall(PetscViewerASCIIPrintf(viewer,"Store in checkpoint number %d (located on disk)\n",rctx->check));
858       break;
859     case 8:
860       PetscCall(PetscViewerASCIIPrintf(viewer,"Restore in checkpoint number %d (located on disk)\n",rctx->check));
861       break;
862     case -1:
863       PetscCall(PetscViewerASCIIPrintf(viewer,"Error!"));
864       break;
865   }
866   PetscFunctionReturn(0);
867 }
868 
869 static PetscErrorCode printwhattodo2(PetscViewer viewer,PetscRevolveInt whattodo,RevolveCTX *rctx,PetscRevolveInt shift)
870 {
871   PetscFunctionBegin;
872   if (!viewer) PetscFunctionReturn(0);
873 
874   switch(whattodo) {
875     case 1:
876       PetscCall(PetscViewerASCIIPrintf(viewer,"[Top Level] Advance from stride %d to stride %d\n",rctx->oldcapo+shift,rctx->capo+shift));
877       break;
878     case 2:
879       PetscCall(PetscViewerASCIIPrintf(viewer,"[Top Level] Store in checkpoint number %d\n",rctx->check));
880       break;
881     case 3:
882       PetscCall(PetscViewerASCIIPrintf(viewer,"[Top Level] First turn: Initialize adjoints and reverse first stride\n"));
883       break;
884     case 4:
885       PetscCall(PetscViewerASCIIPrintf(viewer,"[Top Level] Forward and reverse one stride\n"));
886       break;
887     case 5:
888       PetscCall(PetscViewerASCIIPrintf(viewer,"[Top Level] Restore in checkpoint number %d\n",rctx->check));
889       break;
890     case 7:
891       PetscCall(PetscViewerASCIIPrintf(viewer,"[Top Level] Store in top-level checkpoint number %d\n",rctx->check));
892       break;
893     case 8:
894       PetscCall(PetscViewerASCIIPrintf(viewer,"[Top Level] Restore in top-level checkpoint number %d\n",rctx->check));
895       break;
896     case -1:
897       PetscCall(PetscViewerASCIIPrintf(viewer,"[Top Level] Error!"));
898       break;
899   }
900   PetscFunctionReturn(0);
901 }
902 
903 static PetscErrorCode InitRevolve(PetscInt fine,PetscInt snaps,RevolveCTX *rctx)
904 {
905   PetscRevolveInt rsnaps,rfine;
906 
907   PetscFunctionBegin;
908   PetscCall(PetscRevolveIntCast(snaps,&rsnaps));
909   PetscCall(PetscRevolveIntCast(fine,&rfine));
910   revolve_reset();
911   revolve_create_offline(rfine,rsnaps);
912   rctx->snaps_in       = rsnaps;
913   rctx->fine           = rfine;
914   rctx->check          = 0;
915   rctx->capo           = 0;
916   rctx->reverseonestep = PETSC_FALSE;
917   /* check stepsleft? */
918   PetscFunctionReturn(0);
919 }
920 
921 static PetscErrorCode FastForwardRevolve(RevolveCTX *rctx)
922 {
923   PetscRevolveInt whattodo;
924 
925   PetscFunctionBegin;
926   whattodo = 0;
927   while (whattodo!=3) { /* we have to fast forward revolve to the beginning of the backward sweep due to unfriendly revolve interface */
928     whattodo = revolve_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where);
929   }
930   PetscFunctionReturn(0);
931 }
932 
933 static PetscErrorCode ApplyRevolve(PetscViewer viewer,SchedulerType stype,RevolveCTX *rctx,PetscRevolveInt total_steps,PetscRevolveInt stepnum,PetscRevolveInt localstepnum,PetscBool toplevel,PetscInt *store)
934 {
935   PetscRevolveInt shift,whattodo;
936 
937   PetscFunctionBegin;
938   *store = 0;
939   if (rctx->stepsleft > 0) { /* advance the solution without checkpointing anything as Revolve requires */
940     rctx->stepsleft--;
941     PetscFunctionReturn(0);
942   }
943   /* let Revolve determine what to do next */
944   shift         = stepnum-localstepnum;
945   rctx->oldcapo = rctx->capo;
946   rctx->capo    = localstepnum;
947 
948   if (!toplevel) whattodo = revolve_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where);
949   else whattodo = revolve2_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where);
950   if (stype == REVOLVE_ONLINE && whattodo == 8) whattodo = 5;
951   if (stype == REVOLVE_ONLINE && whattodo == 7) whattodo = 2;
952   if (!toplevel) PetscCall(printwhattodo(viewer,whattodo,rctx,shift));
953   else PetscCall(printwhattodo2(viewer,whattodo,rctx,shift));
954   PetscCheck(whattodo != -1,PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in the Revolve library");
955   if (whattodo == 1) { /* advance some time steps */
956     if (stype == REVOLVE_ONLINE && rctx->capo >= total_steps-1) {
957       revolve_turn(total_steps,&rctx->capo,&rctx->fine);
958       if (!toplevel) PetscCall(printwhattodo(viewer,whattodo,rctx,shift));
959       else PetscCall(printwhattodo2(viewer,whattodo,rctx,shift));
960     }
961     rctx->stepsleft = rctx->capo-rctx->oldcapo-1;
962   }
963   if (whattodo == 3 || whattodo == 4) { /* ready for a reverse step */
964     rctx->reverseonestep = PETSC_TRUE;
965   }
966   if (whattodo == 5) { /* restore a checkpoint and ask Revolve what to do next */
967     rctx->oldcapo = rctx->capo;
968     if (!toplevel) whattodo = revolve_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where); /* must return 1 or 3 or 4*/
969     else whattodo = revolve2_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where);
970     if (!toplevel) PetscCall(printwhattodo(viewer,whattodo,rctx,shift));
971     else PetscCall(printwhattodo2(viewer,whattodo,rctx,shift));
972     if (whattodo == 3 || whattodo == 4) rctx->reverseonestep = PETSC_TRUE;
973     if (whattodo == 1) rctx->stepsleft = rctx->capo-rctx->oldcapo;
974   }
975   if (whattodo == 7) { /* save the checkpoint to disk */
976     *store = 2;
977     rctx->oldcapo = rctx->capo;
978     whattodo = revolve_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where); /* must return 1 */
979     PetscCall(printwhattodo(viewer,whattodo,rctx,shift));
980     rctx->stepsleft = rctx->capo-rctx->oldcapo-1;
981   }
982   if (whattodo == 2) { /* store a checkpoint to RAM and ask Revolve how many time steps to advance next */
983     *store = 1;
984     rctx->oldcapo = rctx->capo;
985     if (!toplevel) whattodo = revolve_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where); /* must return 1 */
986     else whattodo = revolve2_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where);
987     if (!toplevel) PetscCall(printwhattodo(viewer,whattodo,rctx,shift));
988     else PetscCall(printwhattodo2(viewer,whattodo,rctx,shift));
989     if (stype == REVOLVE_ONLINE && rctx->capo >= total_steps-1) {
990       revolve_turn(total_steps,&rctx->capo,&rctx->fine);
991       PetscCall(printwhattodo(viewer,whattodo,rctx,shift));
992     }
993     rctx->stepsleft = rctx->capo-rctx->oldcapo-1;
994   }
995   PetscFunctionReturn(0);
996 }
997 
998 static PetscErrorCode TSTrajectoryMemorySet_ROF(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X)
999 {
1000   Stack           *stack = &tjsch->stack;
1001   PetscInt        store;
1002   StackElement    e;
1003   PetscRevolveInt rtotal_steps,rstepnum;
1004   CheckpointType  cptype;
1005 
1006   PetscFunctionBegin;
1007   if (!stack->solution_only && stepnum == 0) PetscFunctionReturn(0);
1008   if (stack->solution_only && stepnum == tjsch->total_steps) PetscFunctionReturn(0);
1009   PetscCall(PetscRevolveIntCast(tjsch->total_steps,&rtotal_steps));
1010   PetscCall(PetscRevolveIntCast(stepnum,&rstepnum));
1011   PetscCall(ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,rtotal_steps,rstepnum,rstepnum,PETSC_FALSE,&store));
1012   if (store == 1) {
1013     PetscCheck(stepnum >= stack->top,PetscObjectComm((PetscObject)ts),PETSC_ERR_MEMC,"Illegal modification of a non-top stack element");
1014     cptype = stack->solution_only ? SOLUTIONONLY : SOLUTION_STAGES;
1015     PetscCall(ElementCreate(ts,cptype,stack,&e));
1016     PetscCall(ElementSet(ts,stack,&e,stepnum,time,X));
1017     PetscCall(StackPush(stack,e));
1018   }
1019   PetscFunctionReturn(0);
1020 }
1021 
1022 static PetscErrorCode TSTrajectoryMemoryGet_ROF(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum)
1023 {
1024   Stack           *stack = &tjsch->stack;
1025   PetscInt        store;
1026   PetscRevolveInt whattodo,shift,rtotal_steps,rstepnum;
1027   StackElement    e;
1028 
1029   PetscFunctionBegin;
1030   if (stepnum == 0 || stepnum == tjsch->total_steps) {
1031     PetscCall(TurnBackward(ts));
1032     tjsch->rctx->reverseonestep = PETSC_FALSE;
1033     PetscFunctionReturn(0);
1034   }
1035   /* restore a checkpoint */
1036   PetscCall(StackTop(stack,&e));
1037   PetscCall(UpdateTS(ts,stack,e,stepnum,PETSC_TRUE));
1038   PetscCall(PetscRevolveIntCast(tjsch->total_steps,&rtotal_steps));
1039   PetscCall(PetscRevolveIntCast(stepnum,&rstepnum));
1040   if (stack->solution_only) { /* start with restoring a checkpoint */
1041     tjsch->rctx->capo = rstepnum;
1042     tjsch->rctx->oldcapo = tjsch->rctx->capo;
1043     shift = 0;
1044     whattodo = revolve_action(&tjsch->rctx->check,&tjsch->rctx->capo,&tjsch->rctx->fine,tjsch->rctx->snaps_in,&tjsch->rctx->info,&tjsch->rctx->where);
1045     PetscCall(printwhattodo(tj->monitor,whattodo,tjsch->rctx,shift));
1046   } else { /* 2 revolve actions: restore a checkpoint and then advance */
1047     PetscCall(ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,rtotal_steps,rstepnum,rstepnum,PETSC_FALSE,&store));
1048     if (tj->monitor) {
1049       PetscCall(PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel));
1050       PetscCall(PetscViewerASCIIPrintf(tj->monitor,"Skip the step from %d to %d (stage values already checkpointed)\n",tjsch->rctx->oldcapo,tjsch->rctx->oldcapo+1));
1051       PetscCall(PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel));
1052     }
1053     if (!tjsch->rctx->reverseonestep && tjsch->rctx->stepsleft > 0) tjsch->rctx->stepsleft--;
1054   }
1055   if (stack->solution_only || (!stack->solution_only && e->stepnum < stepnum)) {
1056     PetscCall(TurnForward(ts));
1057     PetscCall(ReCompute(ts,tjsch,e->stepnum,stepnum));
1058   }
1059   if ((stack->solution_only && e->stepnum+1 == stepnum) || (!stack->solution_only && e->stepnum == stepnum)) {
1060     PetscCall(StackPop(stack,&e));
1061   }
1062   tjsch->rctx->reverseonestep = PETSC_FALSE;
1063   PetscFunctionReturn(0);
1064 }
1065 
1066 static PetscErrorCode TSTrajectoryMemorySet_RON(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X)
1067 {
1068   Stack           *stack = &tjsch->stack;
1069   Vec             *Y;
1070   PetscInt        i,store;
1071   PetscReal       timeprev;
1072   StackElement    e;
1073   RevolveCTX      *rctx = tjsch->rctx;
1074   PetscRevolveInt rtotal_steps,rstepnum;
1075   CheckpointType  cptype;
1076 
1077   PetscFunctionBegin;
1078   if (!stack->solution_only && stepnum == 0) PetscFunctionReturn(0);
1079   if (stack->solution_only && stepnum == tjsch->total_steps) PetscFunctionReturn(0);
1080   PetscCall(PetscRevolveIntCast(tjsch->total_steps,&rtotal_steps));
1081   PetscCall(PetscRevolveIntCast(stepnum,&rstepnum));
1082   PetscCall(ApplyRevolve(tj->monitor,tjsch->stype,rctx,rtotal_steps,rstepnum,rstepnum,PETSC_FALSE,&store));
1083   if (store == 1) {
1084     if (rctx->check != stack->top+1) { /* overwrite some non-top checkpoint in the stack */
1085       PetscCall(StackFind(stack,&e,rctx->check));
1086       if (HaveSolution(e->cptype)) {
1087         PetscCall(VecCopy(X,e->X));
1088       }
1089       if (HaveStages(e->cptype)) {
1090         PetscCall(TSGetStages(ts,&stack->numY,&Y));
1091         for (i=0;i<stack->numY;i++) {
1092           PetscCall(VecCopy(Y[i],e->Y[i]));
1093         }
1094       }
1095       e->stepnum  = stepnum;
1096       e->time     = time;
1097       PetscCall(TSGetPrevTime(ts,&timeprev));
1098       e->timeprev = timeprev;
1099     } else {
1100       PetscCheck(stepnum >= stack->top,PetscObjectComm((PetscObject)ts),PETSC_ERR_MEMC,"Illegal modification of a non-top stack element");
1101       cptype = stack->solution_only ? SOLUTIONONLY : SOLUTION_STAGES;
1102       PetscCall(ElementCreate(ts,cptype,stack,&e));
1103       PetscCall(ElementSet(ts,stack,&e,stepnum,time,X));
1104       PetscCall(StackPush(stack,e));
1105     }
1106   }
1107   PetscFunctionReturn(0);
1108 }
1109 
1110 static PetscErrorCode TSTrajectoryMemoryGet_RON(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum)
1111 {
1112   Stack           *stack = &tjsch->stack;
1113   PetscRevolveInt whattodo,shift,rstepnum;
1114   StackElement    e;
1115 
1116   PetscFunctionBegin;
1117   if (stepnum == 0 || stepnum == tjsch->total_steps) {
1118     PetscCall(TurnBackward(ts));
1119     tjsch->rctx->reverseonestep = PETSC_FALSE;
1120     PetscFunctionReturn(0);
1121   }
1122   PetscCall(PetscRevolveIntCast(stepnum,&rstepnum));
1123   tjsch->rctx->capo = rstepnum;
1124   tjsch->rctx->oldcapo = tjsch->rctx->capo;
1125   shift = 0;
1126   whattodo = revolve_action(&tjsch->rctx->check,&tjsch->rctx->capo,&tjsch->rctx->fine,tjsch->rctx->snaps_in,&tjsch->rctx->info,&tjsch->rctx->where); /* whattodo=restore */
1127   if (whattodo == 8) whattodo = 5;
1128   PetscCall(printwhattodo(tj->monitor,whattodo,tjsch->rctx,shift));
1129   /* restore a checkpoint */
1130   PetscCall(StackFind(stack,&e,tjsch->rctx->check));
1131   PetscCall(UpdateTS(ts,stack,e,stepnum,PETSC_TRUE));
1132   if (!stack->solution_only) { /* whattodo must be 5 */
1133     /* ask Revolve what to do next */
1134     tjsch->rctx->oldcapo = tjsch->rctx->capo;
1135     whattodo = revolve_action(&tjsch->rctx->check,&tjsch->rctx->capo,&tjsch->rctx->fine,tjsch->rctx->snaps_in,&tjsch->rctx->info,&tjsch->rctx->where); /* must return 1 or 3 or 4*/
1136     PetscCall(printwhattodo(tj->monitor,whattodo,tjsch->rctx,shift));
1137     if (whattodo == 3 || whattodo == 4) tjsch->rctx->reverseonestep = PETSC_TRUE;
1138     if (whattodo == 1) tjsch->rctx->stepsleft = tjsch->rctx->capo-tjsch->rctx->oldcapo;
1139     if (tj->monitor) {
1140       PetscCall(PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel));
1141       PetscCall(PetscViewerASCIIPrintf(tj->monitor,"Skip the step from %d to %d (stage values already checkpointed)\n",tjsch->rctx->oldcapo,tjsch->rctx->oldcapo+1));
1142       PetscCall(PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel));
1143     }
1144     if (!tjsch->rctx->reverseonestep && tjsch->rctx->stepsleft > 0) tjsch->rctx->stepsleft--;
1145   }
1146   if (stack->solution_only || (!stack->solution_only && e->stepnum < stepnum)) {
1147     PetscCall(TurnForward(ts));
1148     PetscCall(ReCompute(ts,tjsch,e->stepnum,stepnum));
1149   }
1150   tjsch->rctx->reverseonestep = PETSC_FALSE;
1151   PetscFunctionReturn(0);
1152 }
1153 
1154 static PetscErrorCode TSTrajectoryMemorySet_TLR(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X)
1155 {
1156   Stack           *stack = &tjsch->stack;
1157   PetscInt        store,localstepnum,laststridesize;
1158   StackElement    e;
1159   PetscBool       done = PETSC_FALSE;
1160   PetscRevolveInt rtotal_steps,rstepnum,rlocalstepnum;
1161   CheckpointType  cptype;
1162 
1163   PetscFunctionBegin;
1164   if (!stack->solution_only && stepnum == 0) PetscFunctionReturn(0);
1165   if (stack->solution_only && stepnum == tjsch->total_steps) PetscFunctionReturn(0);
1166 
1167   localstepnum = stepnum%tjsch->stride;
1168   laststridesize = tjsch->total_steps%tjsch->stride;
1169   if (!laststridesize) laststridesize = tjsch->stride;
1170 
1171   if (!tjsch->recompute) {
1172     PetscCall(TopLevelStore(tj,ts,tjsch,stepnum,localstepnum,laststridesize,&done));
1173     /* revolve is needed for the last stride; different starting points for last stride between solutin_only and !solutin_only */
1174     if (!stack->solution_only && !tjsch->save_stack && stepnum <= tjsch->total_steps-laststridesize) PetscFunctionReturn(0);
1175     if (stack->solution_only && !tjsch->save_stack && stepnum < tjsch->total_steps-laststridesize) PetscFunctionReturn(0);
1176   }
1177   if (tjsch->save_stack && done) {
1178     PetscCall(InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx));
1179     PetscFunctionReturn(0);
1180   }
1181   if (laststridesize < tjsch->stride) {
1182     if (stack->solution_only && stepnum == tjsch->total_steps-laststridesize && !tjsch->recompute) { /* step tjsch->total_steps-laststridesize-1 is skipped, but the next step is not */
1183       PetscCall(InitRevolve(laststridesize,tjsch->max_cps_ram,tjsch->rctx));
1184     }
1185     if (!stack->solution_only && stepnum == tjsch->total_steps-laststridesize+1 && !tjsch->recompute) { /* step tjsch->total_steps-laststridesize is skipped, but the next step is not */
1186       PetscCall(InitRevolve(laststridesize,tjsch->max_cps_ram,tjsch->rctx));
1187     }
1188   }
1189   PetscCall(PetscRevolveIntCast(tjsch->total_steps,&rtotal_steps));
1190   PetscCall(PetscRevolveIntCast(stepnum,&rstepnum));
1191   PetscCall(PetscRevolveIntCast(localstepnum,&rlocalstepnum));
1192   PetscCall(ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,rtotal_steps,rstepnum,rlocalstepnum,PETSC_FALSE,&store));
1193   if (store == 1) {
1194     PetscCheck(localstepnum >= stack->top,PetscObjectComm((PetscObject)ts),PETSC_ERR_MEMC,"Illegal modification of a non-top stack element");
1195     cptype = stack->solution_only ? SOLUTIONONLY : SOLUTION_STAGES;
1196     PetscCall(ElementCreate(ts,cptype,stack,&e));
1197     PetscCall(ElementSet(ts,stack,&e,stepnum,time,X));
1198     PetscCall(StackPush(stack,e));
1199   }
1200   PetscFunctionReturn(0);
1201 }
1202 
1203 static PetscErrorCode TSTrajectoryMemoryGet_TLR(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum)
1204 {
1205   Stack           *stack = &tjsch->stack;
1206   PetscRevolveInt whattodo,shift,rstepnum,rlocalstepnum,rtotal_steps;
1207   PetscInt        localstepnum,stridenum,laststridesize,store;
1208   StackElement    e;
1209   CheckpointType  cptype;
1210 
1211   PetscFunctionBegin;
1212   localstepnum = stepnum%tjsch->stride;
1213   stridenum    = stepnum/tjsch->stride;
1214   if (stepnum == tjsch->total_steps) {
1215     PetscCall(TurnBackward(ts));
1216     tjsch->rctx->reverseonestep = PETSC_FALSE;
1217     PetscFunctionReturn(0);
1218   }
1219   laststridesize = tjsch->total_steps%tjsch->stride;
1220   if (!laststridesize) laststridesize = tjsch->stride;
1221   PetscCall(PetscRevolveIntCast(tjsch->total_steps,&rtotal_steps));
1222   PetscCall(PetscRevolveIntCast(stepnum,&rstepnum));
1223   PetscCall(PetscRevolveIntCast(localstepnum,&rlocalstepnum));
1224   if (stack->solution_only) {
1225     /* fill stack */
1226     if (localstepnum == 0 && stepnum <= tjsch->total_steps-laststridesize) {
1227       if (tjsch->save_stack) {
1228         PetscCall(StackLoadAll(tj,ts,stack,stridenum));
1229         PetscCall(InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx));
1230         PetscCall(FastForwardRevolve(tjsch->rctx));
1231         tjsch->skip_trajectory = PETSC_TRUE;
1232         PetscCall(TurnForward(ts));
1233         PetscCall(ReCompute(ts,tjsch,stridenum*tjsch->stride-1,stridenum*tjsch->stride));
1234         tjsch->skip_trajectory = PETSC_FALSE;
1235       } else {
1236         PetscCall(LoadSingle(tj,ts,stack,stridenum));
1237         PetscCall(InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx));
1238         PetscCall(TurnForward(ts));
1239         PetscCall(ReCompute(ts,tjsch,(stridenum-1)*tjsch->stride,stridenum*tjsch->stride));
1240       }
1241       PetscFunctionReturn(0);
1242     }
1243     /* restore a checkpoint */
1244     PetscCall(StackTop(stack,&e));
1245     PetscCall(UpdateTS(ts,stack,e,stepnum,PETSC_TRUE));
1246     /* start with restoring a checkpoint */
1247     tjsch->rctx->capo = rstepnum;
1248     tjsch->rctx->oldcapo = tjsch->rctx->capo;
1249     shift = rstepnum-rlocalstepnum;
1250     whattodo = revolve_action(&tjsch->rctx->check,&tjsch->rctx->capo,&tjsch->rctx->fine,tjsch->rctx->snaps_in,&tjsch->rctx->info,&tjsch->rctx->where);
1251     PetscCall(printwhattodo(tj->monitor,whattodo,tjsch->rctx,shift));
1252     PetscCall(TurnForward(ts));
1253     PetscCall(ReCompute(ts,tjsch,e->stepnum,stepnum));
1254     if (e->stepnum+1 == stepnum) {
1255       PetscCall(StackPop(stack,&e));
1256     }
1257   } else {
1258     /* fill stack with info */
1259     if (localstepnum == 0 && tjsch->total_steps-stepnum >= laststridesize) {
1260       if (tjsch->save_stack) {
1261         PetscCall(StackLoadAll(tj,ts,stack,stridenum));
1262         PetscCall(InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx));
1263         PetscCall(FastForwardRevolve(tjsch->rctx));
1264       } else {
1265         PetscRevolveInt rnum;
1266         PetscCall(LoadSingle(tj,ts,stack,stridenum));
1267         PetscCall(InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx));
1268         PetscCall(PetscRevolveIntCast((stridenum-1)*tjsch->stride+1,&rnum));
1269         PetscCall(ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,rtotal_steps,rnum,1,PETSC_FALSE,&store));
1270         if (tj->monitor) {
1271           PetscCall(PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel));
1272           PetscCall(PetscViewerASCIIPrintf(tj->monitor,"Skip the step from %" PetscInt_FMT " to %" PetscInt_FMT " (stage values already checkpointed)\n",(stridenum-1)*tjsch->stride+tjsch->rctx->oldcapo,(stridenum-1)*tjsch->stride+tjsch->rctx->oldcapo+1));
1273           PetscCall(PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel));
1274         }
1275         cptype = SOLUTION_STAGES;
1276         PetscCall(ElementCreate(ts,cptype,stack,&e));
1277         PetscCall(ElementSet(ts,stack,&e,(stridenum-1)*tjsch->stride+1,ts->ptime,ts->vec_sol));
1278         PetscCall(StackPush(stack,e));
1279         PetscCall(TurnForward(ts));
1280         PetscCall(ReCompute(ts,tjsch,e->stepnum,stridenum*tjsch->stride));
1281       }
1282       PetscFunctionReturn(0);
1283     }
1284     /* restore a checkpoint */
1285     PetscCall(StackTop(stack,&e));
1286     PetscCall(UpdateTS(ts,stack,e,stepnum,PETSC_TRUE));
1287     /* 2 revolve actions: restore a checkpoint and then advance */
1288     PetscCall(ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,rtotal_steps,rstepnum,rlocalstepnum,PETSC_FALSE,&store));
1289     if (tj->monitor) {
1290       PetscCall(PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel));
1291       PetscCall(PetscViewerASCIIPrintf(tj->monitor,"Skip the step from %" PetscInt_FMT " to %" PetscInt_FMT " (stage values already checkpointed)\n",stepnum-localstepnum+tjsch->rctx->oldcapo,stepnum-localstepnum+tjsch->rctx->oldcapo+1));
1292       PetscCall(PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel));
1293     }
1294     if (!tjsch->rctx->reverseonestep && tjsch->rctx->stepsleft > 0) tjsch->rctx->stepsleft--;
1295     if (e->stepnum < stepnum) {
1296       PetscCall(TurnForward(ts));
1297       PetscCall(ReCompute(ts,tjsch,e->stepnum,stepnum));
1298     }
1299     if (e->stepnum == stepnum) {
1300       PetscCall(StackPop(stack,&e));
1301     }
1302   }
1303   tjsch->rctx->reverseonestep = PETSC_FALSE;
1304   PetscFunctionReturn(0);
1305 }
1306 
1307 static PetscErrorCode TSTrajectoryMemorySet_TLTR(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X)
1308 {
1309   Stack           *stack = &tjsch->stack;
1310   PetscInt        store,localstepnum,stridenum,laststridesize;
1311   StackElement    e;
1312   PetscBool       done = PETSC_FALSE;
1313   PetscRevolveInt rlocalstepnum,rstepnum,rtotal_steps;
1314 
1315   PetscFunctionBegin;
1316   if (!stack->solution_only && stepnum == 0) PetscFunctionReturn(0);
1317   if (stack->solution_only && stepnum == tjsch->total_steps) PetscFunctionReturn(0);
1318 
1319   localstepnum = stepnum%tjsch->stride; /* index at the bottom level (inside a stride) */
1320   stridenum    = stepnum/tjsch->stride; /* index at the top level */
1321   laststridesize = tjsch->total_steps%tjsch->stride;
1322   if (!laststridesize) laststridesize = tjsch->stride;
1323   if (stack->solution_only && localstepnum == 0 && !tjsch->rctx2->reverseonestep) {
1324     PetscCall(PetscRevolveIntCast((tjsch->total_steps+tjsch->stride-1)/tjsch->stride,&rtotal_steps));
1325     PetscCall(PetscRevolveIntCast(stridenum,&rstepnum));
1326     PetscCall(ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx2,rtotal_steps,rstepnum,rstepnum,PETSC_TRUE,&tjsch->store_stride));
1327     if (laststridesize < tjsch->stride && stepnum == tjsch->total_steps-laststridesize) {
1328       PetscCall(InitRevolve(laststridesize,tjsch->max_cps_ram,tjsch->rctx));
1329     }
1330   }
1331   if (!stack->solution_only && localstepnum == 1 && !tjsch->rctx2->reverseonestep) {
1332     PetscCall(PetscRevolveIntCast((tjsch->total_steps+tjsch->stride-1)/tjsch->stride,&rtotal_steps));
1333     PetscCall(PetscRevolveIntCast(stridenum,&rstepnum));
1334     PetscCall(ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx2,rtotal_steps,rstepnum,rstepnum,PETSC_TRUE,&tjsch->store_stride));
1335     if (laststridesize < tjsch->stride && stepnum == tjsch->total_steps-laststridesize+1) {
1336       PetscCall(InitRevolve(laststridesize,tjsch->max_cps_ram,tjsch->rctx));
1337     }
1338   }
1339   if (tjsch->store_stride) {
1340     PetscCall(TopLevelStore(tj,ts,tjsch,stepnum,localstepnum,laststridesize,&done));
1341     if (done) {
1342       PetscCall(InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx));
1343       PetscFunctionReturn(0);
1344     }
1345   }
1346   if (stepnum < tjsch->total_steps-laststridesize) {
1347     if (tjsch->save_stack && !tjsch->store_stride && !tjsch->rctx2->reverseonestep) PetscFunctionReturn(0); /* store or forward-and-reverse at top level trigger revolve at bottom level */
1348     if (!tjsch->save_stack && !tjsch->rctx2->reverseonestep) PetscFunctionReturn(0); /* store operation does not require revolve be called at bottom level */
1349   }
1350   /* Skipping stepnum=0 for !stack->only is enough for TLR, but not for TLTR. Here we skip the first step for each stride so that the top-level revolve is applied (always at localstepnum=1) ahead of the bottom-level revolve */
1351   if (!stack->solution_only && localstepnum == 0 && stepnum != tjsch->total_steps && !tjsch->recompute) PetscFunctionReturn(0);
1352   PetscCall(PetscRevolveIntCast(tjsch->total_steps,&rtotal_steps));
1353   PetscCall(PetscRevolveIntCast(stepnum,&rstepnum));
1354   PetscCall(PetscRevolveIntCast(localstepnum,&rlocalstepnum));
1355   PetscCall(ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,rtotal_steps,rstepnum,rlocalstepnum,PETSC_FALSE,&store));
1356   if (store == 1) {
1357     CheckpointType cptype;
1358     PetscCheck(localstepnum >= stack->top,PetscObjectComm((PetscObject)ts),PETSC_ERR_MEMC,"Illegal modification of a non-top stack element");
1359     cptype = stack->solution_only ? SOLUTIONONLY : SOLUTION_STAGES;
1360     PetscCall(ElementCreate(ts,cptype,stack,&e));
1361     PetscCall(ElementSet(ts,stack,&e,stepnum,time,X));
1362     PetscCall(StackPush(stack,e));
1363   }
1364   PetscFunctionReturn(0);
1365 }
1366 
1367 static PetscErrorCode TSTrajectoryMemoryGet_TLTR(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum)
1368 {
1369   Stack           *stack = &tjsch->stack;
1370   DiskStack       *diskstack = &tjsch->diskstack;
1371   PetscInt        localstepnum,stridenum,restoredstridenum,laststridesize,store;
1372   StackElement    e;
1373   PetscRevolveInt whattodo,shift;
1374   PetscRevolveInt rtotal_steps,rstepnum,rlocalstepnum;
1375 
1376   PetscFunctionBegin;
1377   localstepnum = stepnum%tjsch->stride;
1378   stridenum    = stepnum/tjsch->stride;
1379   if (stepnum == tjsch->total_steps) {
1380     PetscCall(TurnBackward(ts));
1381     tjsch->rctx->reverseonestep = PETSC_FALSE;
1382     PetscFunctionReturn(0);
1383   }
1384   laststridesize = tjsch->total_steps%tjsch->stride;
1385   if (!laststridesize) laststridesize = tjsch->stride;
1386   /*
1387    Last stride can be adjoined directly. All the other strides require that the stack in memory be ready before an adjoint step is taken (at the end of each stride). The following two cases need to be addressed differently:
1388      Case 1 (save_stack)
1389        Restore a disk checkpoint; update TS with the last element in the restored data; recompute to the current point.
1390      Case 2 (!save_stack)
1391        Restore a disk checkpoint; update TS with the restored point; recompute to the current point.
1392   */
1393   if (localstepnum == 0 && stepnum <= tjsch->total_steps-laststridesize) {
1394     /* restore the top element in the stack for disk checkpoints */
1395     restoredstridenum = diskstack->container[diskstack->top];
1396     tjsch->rctx2->reverseonestep = PETSC_FALSE;
1397     /* top-level revolve must be applied before current step, just like the solution_only mode for single-level revolve */
1398     if (!tjsch->save_stack && stack->solution_only) { /* start with restoring a checkpoint */
1399       PetscCall(PetscRevolveIntCast(stridenum,&rstepnum));
1400       tjsch->rctx2->capo = rstepnum;
1401       tjsch->rctx2->oldcapo = tjsch->rctx2->capo;
1402       shift = 0;
1403       whattodo = revolve2_action(&tjsch->rctx2->check,&tjsch->rctx2->capo,&tjsch->rctx2->fine,tjsch->rctx2->snaps_in,&tjsch->rctx2->info,&tjsch->rctx2->where);
1404       PetscCall(printwhattodo2(tj->monitor,whattodo,tjsch->rctx2,shift));
1405     } else { /* 2 revolve actions: restore a checkpoint and then advance */
1406       PetscCall(PetscRevolveIntCast((tjsch->total_steps+tjsch->stride-1)/tjsch->stride,&rtotal_steps));
1407       PetscCall(PetscRevolveIntCast(stridenum,&rstepnum));
1408       PetscCall(ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx2,rtotal_steps,rstepnum,rstepnum,PETSC_TRUE,&tjsch->store_stride));
1409       if (tj->monitor) {
1410         PetscCall(PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel));
1411         PetscCall(PetscViewerASCIIPrintf(tj->monitor,"[Top Level] Skip the stride from %d to %d (stage values already checkpointed)\n",tjsch->rctx2->oldcapo,tjsch->rctx2->oldcapo+1));
1412         PetscCall(PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel));
1413       }
1414       if (!tjsch->rctx2->reverseonestep && tjsch->rctx2->stepsleft > 0) tjsch->rctx2->stepsleft--;
1415     }
1416     /* fill stack */
1417     if (stack->solution_only) {
1418       if (tjsch->save_stack) {
1419         if (restoredstridenum < stridenum) {
1420           PetscCall(StackLoadLast(tj,ts,stack,restoredstridenum));
1421         } else {
1422           PetscCall(StackLoadAll(tj,ts,stack,restoredstridenum));
1423         }
1424         /* recompute one step ahead */
1425         tjsch->skip_trajectory = PETSC_TRUE;
1426         PetscCall(TurnForward(ts));
1427         PetscCall(ReCompute(ts,tjsch,stridenum*tjsch->stride-1,stridenum*tjsch->stride));
1428         tjsch->skip_trajectory = PETSC_FALSE;
1429         if (restoredstridenum < stridenum) {
1430           PetscCall(InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx));
1431           PetscCall(TurnForward(ts));
1432           PetscCall(ReCompute(ts,tjsch,restoredstridenum*tjsch->stride,stepnum));
1433         } else { /* stack ready, fast forward revolve status */
1434           PetscCall(InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx));
1435           PetscCall(FastForwardRevolve(tjsch->rctx));
1436         }
1437       } else {
1438         PetscCall(LoadSingle(tj,ts,stack,restoredstridenum));
1439         PetscCall(InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx));
1440         PetscCall(TurnForward(ts));
1441         PetscCall(ReCompute(ts,tjsch,(restoredstridenum-1)*tjsch->stride,stepnum));
1442       }
1443     } else {
1444       if (tjsch->save_stack) {
1445         if (restoredstridenum < stridenum) {
1446           PetscCall(StackLoadLast(tj,ts,stack,restoredstridenum));
1447           /* reset revolve */
1448           PetscCall(InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx));
1449           PetscCall(TurnForward(ts));
1450           PetscCall(ReCompute(ts,tjsch,restoredstridenum*tjsch->stride,stepnum));
1451         } else { /* stack ready, fast forward revolve status */
1452           PetscCall(StackLoadAll(tj,ts,stack,restoredstridenum));
1453           PetscCall(InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx));
1454           PetscCall(FastForwardRevolve(tjsch->rctx));
1455         }
1456       } else {
1457         PetscCall(LoadSingle(tj,ts,stack,restoredstridenum));
1458         PetscCall(InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx));
1459         /* push first element to stack */
1460         if (tjsch->store_stride || tjsch->rctx2->reverseonestep) {
1461           CheckpointType cptype = SOLUTION_STAGES;
1462           shift = (restoredstridenum-1)*tjsch->stride-localstepnum;
1463           PetscCall(PetscRevolveIntCast(tjsch->total_steps,&rtotal_steps));
1464           PetscCall(PetscRevolveIntCast((restoredstridenum-1)*tjsch->stride+1,&rstepnum));
1465           PetscCall(ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,rtotal_steps,rstepnum,1,PETSC_FALSE,&store));
1466           if (tj->monitor) {
1467             PetscCall(PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel));
1468             PetscCall(PetscViewerASCIIPrintf(tj->monitor,"Skip the step from %" PetscInt_FMT " to %" PetscInt_FMT " (stage values already checkpointed)\n",(restoredstridenum-1)*tjsch->stride,(restoredstridenum-1)*tjsch->stride+1));
1469             PetscCall(PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel));
1470           }
1471           PetscCall(ElementCreate(ts,cptype,stack,&e));
1472           PetscCall(ElementSet(ts,stack,&e,(restoredstridenum-1)*tjsch->stride+1,ts->ptime,ts->vec_sol));
1473           PetscCall(StackPush(stack,e));
1474         }
1475         PetscCall(TurnForward(ts));
1476         PetscCall(ReCompute(ts,tjsch,(restoredstridenum-1)*tjsch->stride+1,stepnum));
1477       }
1478     }
1479     if (restoredstridenum == stridenum) diskstack->top--;
1480     tjsch->rctx->reverseonestep = PETSC_FALSE;
1481     PetscFunctionReturn(0);
1482   }
1483 
1484   if (stack->solution_only) {
1485     /* restore a checkpoint */
1486     PetscCall(StackTop(stack,&e));
1487     PetscCall(UpdateTS(ts,stack,e,stepnum,PETSC_TRUE));
1488     /* start with restoring a checkpoint */
1489     PetscCall(PetscRevolveIntCast(stepnum,&rstepnum));
1490     PetscCall(PetscRevolveIntCast(localstepnum,&rlocalstepnum));
1491     tjsch->rctx->capo = rstepnum;
1492     tjsch->rctx->oldcapo = tjsch->rctx->capo;
1493     shift = rstepnum-rlocalstepnum;
1494     whattodo = revolve_action(&tjsch->rctx->check,&tjsch->rctx->capo,&tjsch->rctx->fine,tjsch->rctx->snaps_in,&tjsch->rctx->info,&tjsch->rctx->where);
1495     PetscCall(printwhattodo(tj->monitor,whattodo,tjsch->rctx,shift));
1496     PetscCall(TurnForward(ts));
1497     PetscCall(ReCompute(ts,tjsch,e->stepnum,stepnum));
1498     if (e->stepnum+1 == stepnum) {
1499       PetscCall(StackPop(stack,&e));
1500     }
1501   } else {
1502     PetscRevolveInt rlocalstepnum;
1503     /* restore a checkpoint */
1504     PetscCall(StackTop(stack,&e));
1505     PetscCall(UpdateTS(ts,stack,e,stepnum,PETSC_TRUE));
1506     /* 2 revolve actions: restore a checkpoint and then advance */
1507     PetscCall(PetscRevolveIntCast(tjsch->total_steps,&rtotal_steps));
1508     PetscCall(PetscRevolveIntCast(stridenum,&rstepnum));
1509     PetscCall(PetscRevolveIntCast(localstepnum,&rlocalstepnum));
1510     PetscCall(ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,rtotal_steps,rstepnum,rlocalstepnum,PETSC_FALSE,&store));
1511     if (tj->monitor) {
1512       PetscCall(PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel));
1513       PetscCall(PetscViewerASCIIPrintf(tj->monitor,"Skip the step from %" PetscInt_FMT " to %" PetscInt_FMT " (stage values already checkpointed)\n",stepnum-localstepnum+tjsch->rctx->oldcapo,stepnum-localstepnum+tjsch->rctx->oldcapo+1));
1514       PetscCall(PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel));
1515     }
1516     if (!tjsch->rctx->reverseonestep && tjsch->rctx->stepsleft > 0) tjsch->rctx->stepsleft--;
1517     if (e->stepnum < stepnum) {
1518       PetscCall(TurnForward(ts));
1519       PetscCall(ReCompute(ts,tjsch,e->stepnum,stepnum));
1520     }
1521     if (e->stepnum == stepnum) {
1522       PetscCall(StackPop(stack,&e));
1523     }
1524   }
1525   tjsch->rctx->reverseonestep = PETSC_FALSE;
1526   PetscFunctionReturn(0);
1527 }
1528 
1529 static PetscErrorCode TSTrajectoryMemorySet_RMS(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X)
1530 {
1531   Stack           *stack = &tjsch->stack;
1532   PetscInt        store;
1533   StackElement    e;
1534   PetscRevolveInt rtotal_steps,rstepnum;
1535 
1536   PetscFunctionBegin;
1537   if (!stack->solution_only && stepnum == 0) PetscFunctionReturn(0);
1538   if (stack->solution_only && stepnum == tjsch->total_steps) PetscFunctionReturn(0);
1539   PetscCall(PetscRevolveIntCast(tjsch->total_steps,&rtotal_steps));
1540   PetscCall(PetscRevolveIntCast(stepnum,&rstepnum));
1541   PetscCall(ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,rtotal_steps,rstepnum,rstepnum,PETSC_FALSE,&store));
1542   if (store == 1) {
1543     CheckpointType cptype;
1544     PetscCheck(stepnum >= stack->top,PetscObjectComm((PetscObject)ts),PETSC_ERR_MEMC,"Illegal modification of a non-top stack element");
1545     cptype = stack->solution_only ? SOLUTIONONLY : SOLUTION_STAGES;
1546     PetscCall(ElementCreate(ts,cptype,stack,&e));
1547     PetscCall(ElementSet(ts,stack,&e,stepnum,time,X));
1548     PetscCall(StackPush(stack,e));
1549   } else if (store == 2) {
1550     PetscCall(DumpSingle(tj,ts,stack,tjsch->rctx->check+1));
1551   }
1552   PetscFunctionReturn(0);
1553 }
1554 
1555 static PetscErrorCode TSTrajectoryMemoryGet_RMS(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum)
1556 {
1557   Stack           *stack = &tjsch->stack;
1558   PetscRevolveInt whattodo,shift,rstepnum;
1559   PetscInt        restart;
1560   PetscBool       ondisk;
1561   StackElement    e;
1562 
1563   PetscFunctionBegin;
1564   if (stepnum == 0 || stepnum == tjsch->total_steps) {
1565     PetscCall(TurnBackward(ts));
1566     tjsch->rctx->reverseonestep = PETSC_FALSE;
1567     PetscFunctionReturn(0);
1568   }
1569   PetscCall(PetscRevolveIntCast(stepnum,&rstepnum));
1570   tjsch->rctx->capo = rstepnum;
1571   tjsch->rctx->oldcapo = tjsch->rctx->capo;
1572   shift = 0;
1573   whattodo = revolve_action(&tjsch->rctx->check,&tjsch->rctx->capo,&tjsch->rctx->fine,tjsch->rctx->snaps_in,&tjsch->rctx->info,&tjsch->rctx->where); /* whattodo=restore */
1574   PetscCall(printwhattodo(tj->monitor,whattodo,tjsch->rctx,shift));
1575   /* restore a checkpoint */
1576   restart = tjsch->rctx->capo;
1577   if (!tjsch->rctx->where) {
1578     ondisk = PETSC_TRUE;
1579     PetscCall(LoadSingle(tj,ts,stack,tjsch->rctx->check+1));
1580     PetscCall(TurnBackward(ts));
1581   } else {
1582     ondisk = PETSC_FALSE;
1583     PetscCall(StackTop(stack,&e));
1584     PetscCall(UpdateTS(ts,stack,e,stepnum,PETSC_TRUE));
1585   }
1586   if (!stack->solution_only) { /* whattodo must be 5 or 8 */
1587     /* ask Revolve what to do next */
1588     tjsch->rctx->oldcapo = tjsch->rctx->capo;
1589     whattodo = revolve_action(&tjsch->rctx->check,&tjsch->rctx->capo,&tjsch->rctx->fine,tjsch->rctx->snaps_in,&tjsch->rctx->info,&tjsch->rctx->where); /* must return 1 or 3 or 4*/
1590     PetscCall(printwhattodo(tj->monitor,whattodo,tjsch->rctx,shift));
1591     if (whattodo == 3 || whattodo == 4) tjsch->rctx->reverseonestep = PETSC_TRUE;
1592     if (whattodo == 1) tjsch->rctx->stepsleft = tjsch->rctx->capo-tjsch->rctx->oldcapo;
1593     if (tj->monitor) {
1594       PetscCall(PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel));
1595       PetscCall(PetscViewerASCIIPrintf(tj->monitor,"Skip the step from %d to %d (stage values already checkpointed)\n",tjsch->rctx->oldcapo,tjsch->rctx->oldcapo+1));
1596       PetscCall(PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel));
1597     }
1598     if (!tjsch->rctx->reverseonestep && tjsch->rctx->stepsleft > 0) tjsch->rctx->stepsleft--;
1599     restart++; /* skip one step */
1600   }
1601   if (stack->solution_only || (!stack->solution_only && restart < stepnum)) {
1602     PetscCall(TurnForward(ts));
1603     PetscCall(ReCompute(ts,tjsch,restart,stepnum));
1604   }
1605   if (!ondisk && ( (stack->solution_only && e->stepnum+1 == stepnum) || (!stack->solution_only && e->stepnum == stepnum))) {
1606     PetscCall(StackPop(stack,&e));
1607   }
1608   tjsch->rctx->reverseonestep = PETSC_FALSE;
1609   PetscFunctionReturn(0);
1610 }
1611 #endif
1612 
1613 #if defined(PETSC_HAVE_CAMS)
1614 /* Optimal offline adjoint checkpointing for multistage time integration methods */
1615 static PetscErrorCode TSTrajectoryMemorySet_AOF(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X)
1616 {
1617   Stack        *stack = &tjsch->stack;
1618   StackElement  e;
1619 
1620   PetscFunctionBegin;
1621   /* skip if no checkpoint to use. This also avoids an error when num_units_avail=0  */
1622   if (tjsch->actx->nextcheckpointstep == -1) PetscFunctionReturn(0);
1623   if (stepnum == 0) { /* When placing the first checkpoint, no need to change the units available */
1624     if (stack->solution_only) {
1625       PetscCall(offline_ca(tjsch->actx->lastcheckpointstep,tjsch->actx->num_units_avail,tjsch->actx->endstep,&tjsch->actx->nextcheckpointstep));
1626     } else {
1627       /* First two arguments must be -1 when first time calling cams */
1628       PetscCall(offline_cams(tjsch->actx->lastcheckpointstep,tjsch->actx->lastcheckpointtype,tjsch->actx->num_units_avail,tjsch->actx->endstep,tjsch->actx->num_stages,&tjsch->actx->nextcheckpointstep,&tjsch->actx->nextcheckpointtype));
1629     }
1630   }
1631 
1632   if (stack->solution_only && stepnum == tjsch->total_steps) PetscFunctionReturn(0);
1633 
1634   if (tjsch->actx->nextcheckpointstep == stepnum) {
1635     PetscCheck(stepnum >= stack->top,PetscObjectComm((PetscObject)ts),PETSC_ERR_MEMC,"Illegal modification of a non-top stack element");
1636 
1637     if (tjsch->actx->nextcheckpointtype == 2) { /* solution + stage values */
1638       if (tj->monitor) {
1639         PetscCall(PetscViewerASCIIPrintf(tj->monitor,"Store in checkpoint number %" PetscInt_FMT " with stage values and solution (located in RAM)\n",stepnum));
1640       }
1641       PetscCall(ElementCreate(ts,SOLUTION_STAGES,stack,&e));
1642       PetscCall(ElementSet(ts,stack,&e,stepnum,time,X));
1643     }
1644     if (tjsch->actx->nextcheckpointtype == 1) {
1645       if (tj->monitor) {
1646         PetscCall(PetscViewerASCIIPrintf(tj->monitor,"Store in checkpoint number %" PetscInt_FMT " with stage values (located in RAM)\n",stepnum));
1647       }
1648       PetscCall(ElementCreate(ts,STAGESONLY,stack,&e));
1649       PetscCall(ElementSet(ts,stack,&e,stepnum,time,X));
1650     }
1651     if (tjsch->actx->nextcheckpointtype == 0) { /* solution only */
1652       if (tj->monitor) {
1653         PetscCall(PetscViewerASCIIPrintf(tj->monitor,"Store in checkpoint number %" PetscInt_FMT " (located in RAM)\n",stepnum));
1654       }
1655       PetscCall(ElementCreate(ts,SOLUTIONONLY,stack,&e));
1656       PetscCall(ElementSet(ts,stack,&e,stepnum,time,X));
1657     }
1658     PetscCall(StackPush(stack,e));
1659 
1660     tjsch->actx->lastcheckpointstep = stepnum;
1661     if (stack->solution_only) {
1662       PetscCall(offline_ca(tjsch->actx->lastcheckpointstep,tjsch->actx->num_units_avail,tjsch->actx->endstep,&tjsch->actx->nextcheckpointstep));
1663       tjsch->actx->num_units_avail--;
1664     } else {
1665       tjsch->actx->lastcheckpointtype = tjsch->actx->nextcheckpointtype;
1666       PetscCall(offline_cams(tjsch->actx->lastcheckpointstep,tjsch->actx->lastcheckpointtype,tjsch->actx->num_units_avail,tjsch->actx->endstep,tjsch->actx->num_stages,&tjsch->actx->nextcheckpointstep,&tjsch->actx->nextcheckpointtype));
1667       if (tjsch->actx->lastcheckpointtype == 2) tjsch->actx->num_units_avail -= tjsch->actx->num_stages+1;
1668       if (tjsch->actx->lastcheckpointtype == 1) tjsch->actx->num_units_avail -= tjsch->actx->num_stages;
1669       if (tjsch->actx->lastcheckpointtype == 0) tjsch->actx->num_units_avail--;
1670     }
1671   }
1672   PetscFunctionReturn(0);
1673 }
1674 
1675 static PetscErrorCode TSTrajectoryMemoryGet_AOF(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum)
1676 {
1677   Stack        *stack = &tjsch->stack;
1678   StackElement  e;
1679   PetscInt      estepnum;
1680 
1681   PetscFunctionBegin;
1682   if (stepnum == 0 || stepnum == tjsch->total_steps) {
1683     PetscCall(TurnBackward(ts));
1684     PetscFunctionReturn(0);
1685   }
1686   /* Restore a checkpoint */
1687   PetscCall(StackTop(stack,&e));
1688   estepnum = e->stepnum;
1689   if (estepnum == stepnum && e->cptype == SOLUTIONONLY) { /* discard the checkpoint if not useful (corner case) */
1690     PetscCall(StackPop(stack,&e));
1691     tjsch->actx->num_units_avail++;
1692     PetscCall(StackTop(stack,&e));
1693     estepnum = e->stepnum;
1694   }
1695   /* Update TS with stage values if an adjoint step can be taken immediately */
1696   if (HaveStages(e->cptype)) {
1697     if (tj->monitor) {
1698       PetscCall(PetscViewerASCIIPrintf(tj->monitor,"Restore in checkpoint number %" PetscInt_FMT " with stage values (located in RAM)\n",e->stepnum));
1699     }
1700     if (e->cptype == STAGESONLY) tjsch->actx->num_units_avail += tjsch->actx->num_stages;
1701     if (e->cptype == SOLUTION_STAGES) tjsch->actx->num_units_avail += tjsch->actx->num_stages+1;
1702   } else {
1703     if (tj->monitor) {
1704       PetscCall(PetscViewerASCIIPrintf(tj->monitor,"Restore in checkpoint number %" PetscInt_FMT " (located in RAM)\n",e->stepnum));
1705     }
1706     tjsch->actx->num_units_avail++;
1707   }
1708   PetscCall(UpdateTS(ts,stack,e,stepnum,PETSC_TRUE));
1709   /* Query the scheduler */
1710   tjsch->actx->lastcheckpointstep = estepnum;
1711   tjsch->actx->endstep = stepnum;
1712   if (stack->solution_only) { /* start with restoring a checkpoint */
1713     PetscCall(offline_ca(tjsch->actx->lastcheckpointstep,tjsch->actx->num_units_avail,tjsch->actx->endstep,&tjsch->actx->nextcheckpointstep));
1714   } else { /* 2 revolve actions: restore a checkpoint and then advance */
1715     tjsch->actx->lastcheckpointtype = e->cptype;
1716     PetscCall(offline_cams(tjsch->actx->lastcheckpointstep,tjsch->actx->lastcheckpointtype,tjsch->actx->num_units_avail,tjsch->actx->endstep,tjsch->actx->num_stages,&tjsch->actx->nextcheckpointstep, &tjsch->actx->nextcheckpointtype));
1717   }
1718   /* Discard the checkpoint if not needed, decrease the number of available checkpoints if it still stays in stack */
1719   if (HaveStages(e->cptype)) {
1720     if (estepnum == stepnum) {
1721       PetscCall(StackPop(stack,&e));
1722     } else {
1723       if (e->cptype == STAGESONLY) tjsch->actx->num_units_avail -= tjsch->actx->num_stages;
1724       if (e->cptype == SOLUTION_STAGES) tjsch->actx->num_units_avail -= tjsch->actx->num_stages+1;
1725     }
1726   } else {
1727     if (estepnum+1 == stepnum) {
1728       PetscCall(StackPop(stack,&e));
1729     } else {
1730       tjsch->actx->num_units_avail--;
1731     }
1732   }
1733   /* Recompute from the restored checkpoint */
1734   if (stack->solution_only || (!stack->solution_only && estepnum < stepnum)) {
1735     PetscCall(TurnForward(ts));
1736     PetscCall(ReCompute(ts,tjsch,estepnum,stepnum));
1737   }
1738   PetscFunctionReturn(0);
1739 }
1740 #endif
1741 
1742 static PetscErrorCode TSTrajectorySet_Memory(TSTrajectory tj,TS ts,PetscInt stepnum,PetscReal time,Vec X)
1743 {
1744   TJScheduler *tjsch = (TJScheduler*)tj->data;
1745 
1746   PetscFunctionBegin;
1747   if (!tjsch->recompute) { /* use global stepnum in the forward sweep */
1748     PetscCall(TSGetStepNumber(ts,&stepnum));
1749   }
1750   /* for consistency */
1751   if (!tjsch->recompute && stepnum == 0) ts->ptime_prev = ts->ptime-ts->time_step;
1752   switch (tjsch->stype) {
1753     case NONE:
1754       if (tj->adjoint_solve_mode) {
1755         PetscCall(TSTrajectoryMemorySet_N(ts,tjsch,stepnum,time,X));
1756       } else {
1757         PetscCall(TSTrajectoryMemorySet_N_2(ts,tjsch,stepnum,time,X));
1758       }
1759       break;
1760     case TWO_LEVEL_NOREVOLVE:
1761       PetscCheck(tj->adjoint_solve_mode,PetscObjectComm((PetscObject)tj),PETSC_ERR_SUP,"Not implemented");
1762       PetscCall(TSTrajectoryMemorySet_TLNR(tj,ts,tjsch,stepnum,time,X));
1763       break;
1764 #if defined(PETSC_HAVE_REVOLVE)
1765     case TWO_LEVEL_REVOLVE:
1766       PetscCheck(tj->adjoint_solve_mode,PetscObjectComm((PetscObject)tj),PETSC_ERR_SUP,"Not implemented");
1767       PetscCall(TSTrajectoryMemorySet_TLR(tj,ts,tjsch,stepnum,time,X));
1768       break;
1769     case TWO_LEVEL_TWO_REVOLVE:
1770       PetscCheck(tj->adjoint_solve_mode,PetscObjectComm((PetscObject)tj),PETSC_ERR_SUP,"Not implemented");
1771       PetscCall(TSTrajectoryMemorySet_TLTR(tj,ts,tjsch,stepnum,time,X));
1772       break;
1773     case REVOLVE_OFFLINE:
1774       PetscCheck(tj->adjoint_solve_mode,PetscObjectComm((PetscObject)tj),PETSC_ERR_SUP,"Not implemented");
1775       PetscCall(TSTrajectoryMemorySet_ROF(tj,ts,tjsch,stepnum,time,X));
1776       break;
1777     case REVOLVE_ONLINE:
1778       PetscCheck(tj->adjoint_solve_mode,PetscObjectComm((PetscObject)tj),PETSC_ERR_SUP,"Not implemented");
1779       PetscCall(TSTrajectoryMemorySet_RON(tj,ts,tjsch,stepnum,time,X));
1780       break;
1781     case REVOLVE_MULTISTAGE:
1782       PetscCheck(tj->adjoint_solve_mode,PetscObjectComm((PetscObject)tj),PETSC_ERR_SUP,"Not implemented");
1783       PetscCall(TSTrajectoryMemorySet_RMS(tj,ts,tjsch,stepnum,time,X));
1784       break;
1785 #endif
1786 #if defined(PETSC_HAVE_CAMS)
1787     case CAMS_OFFLINE:
1788       PetscCheck(tj->adjoint_solve_mode,PetscObjectComm((PetscObject)tj),PETSC_ERR_SUP,"Not implemented");
1789       PetscCall(TSTrajectoryMemorySet_AOF(tj,ts,tjsch,stepnum,time,X));
1790       break;
1791 #endif
1792     default:
1793       break;
1794   }
1795   PetscFunctionReturn(0);
1796 }
1797 
1798 static PetscErrorCode TSTrajectoryGet_Memory(TSTrajectory tj,TS ts,PetscInt stepnum,PetscReal *t)
1799 {
1800   TJScheduler *tjsch = (TJScheduler*)tj->data;
1801 
1802   PetscFunctionBegin;
1803   if (tj->adjoint_solve_mode && stepnum == 0) {
1804     PetscCall(TSTrajectoryReset(tj)); /* reset TSTrajectory so users do not need to reset TSTrajectory */
1805     PetscFunctionReturn(0);
1806   }
1807   switch (tjsch->stype) {
1808     case NONE:
1809       if (tj->adjoint_solve_mode) {
1810         PetscCall(TSTrajectoryMemoryGet_N(ts,tjsch,stepnum));
1811       } else {
1812         PetscCall(TSTrajectoryMemoryGet_N_2(ts,tjsch,stepnum));
1813       }
1814       break;
1815     case TWO_LEVEL_NOREVOLVE:
1816       PetscCheck(tj->adjoint_solve_mode,PetscObjectComm((PetscObject)tj),PETSC_ERR_SUP,"Not implemented");
1817       PetscCall(TSTrajectoryMemoryGet_TLNR(tj,ts,tjsch,stepnum));
1818       break;
1819 #if defined(PETSC_HAVE_REVOLVE)
1820     case TWO_LEVEL_REVOLVE:
1821       PetscCheck(tj->adjoint_solve_mode,PetscObjectComm((PetscObject)tj),PETSC_ERR_SUP,"Not implemented");
1822       PetscCall(TSTrajectoryMemoryGet_TLR(tj,ts,tjsch,stepnum));
1823       break;
1824     case TWO_LEVEL_TWO_REVOLVE:
1825       PetscCheck(tj->adjoint_solve_mode,PetscObjectComm((PetscObject)tj),PETSC_ERR_SUP,"Not implemented");
1826       PetscCall(TSTrajectoryMemoryGet_TLTR(tj,ts,tjsch,stepnum));
1827       break;
1828     case REVOLVE_OFFLINE:
1829       PetscCheck(tj->adjoint_solve_mode,PetscObjectComm((PetscObject)tj),PETSC_ERR_SUP,"Not implemented");
1830       PetscCall(TSTrajectoryMemoryGet_ROF(tj,ts,tjsch,stepnum));
1831       break;
1832     case REVOLVE_ONLINE:
1833       PetscCheck(tj->adjoint_solve_mode,PetscObjectComm((PetscObject)tj),PETSC_ERR_SUP,"Not implemented");
1834       PetscCall(TSTrajectoryMemoryGet_RON(tj,ts,tjsch,stepnum));
1835       break;
1836     case REVOLVE_MULTISTAGE:
1837       PetscCheck(tj->adjoint_solve_mode,PetscObjectComm((PetscObject)tj),PETSC_ERR_SUP,"Not implemented");
1838       PetscCall(TSTrajectoryMemoryGet_RMS(tj,ts,tjsch,stepnum));
1839       break;
1840 #endif
1841 #if defined(PETSC_HAVE_CAMS)
1842     case CAMS_OFFLINE:
1843       PetscCheck(tj->adjoint_solve_mode,PetscObjectComm((PetscObject)tj),PETSC_ERR_SUP,"Not implemented");
1844       PetscCall(TSTrajectoryMemoryGet_AOF(tj,ts,tjsch,stepnum));
1845       break;
1846 #endif
1847     default:
1848       break;
1849   }
1850   PetscFunctionReturn(0);
1851 }
1852 
1853 PETSC_UNUSED static PetscErrorCode TSTrajectorySetStride_Memory(TSTrajectory tj,PetscInt stride)
1854 {
1855   TJScheduler *tjsch = (TJScheduler*)tj->data;
1856 
1857   PetscFunctionBegin;
1858   tjsch->stride = stride;
1859   PetscFunctionReturn(0);
1860 }
1861 
1862 static PetscErrorCode TSTrajectorySetMaxCpsRAM_Memory(TSTrajectory tj,PetscInt max_cps_ram)
1863 {
1864   TJScheduler *tjsch = (TJScheduler*)tj->data;
1865 
1866   PetscFunctionBegin;
1867   tjsch->max_cps_ram = max_cps_ram;
1868   PetscFunctionReturn(0);
1869 }
1870 
1871 static PetscErrorCode TSTrajectorySetMaxCpsDisk_Memory(TSTrajectory tj,PetscInt max_cps_disk)
1872 {
1873   TJScheduler *tjsch = (TJScheduler*)tj->data;
1874 
1875   PetscFunctionBegin;
1876   tjsch->max_cps_disk = max_cps_disk;
1877   PetscFunctionReturn(0);
1878 }
1879 
1880 static PetscErrorCode TSTrajectorySetMaxUnitsRAM_Memory(TSTrajectory tj,PetscInt max_units_ram)
1881 {
1882   TJScheduler *tjsch = (TJScheduler*)tj->data;
1883 
1884   PetscFunctionBegin;
1885   PetscCheck(tjsch->max_cps_ram,PetscObjectComm((PetscObject)tj),PETSC_ERR_ARG_INCOMP,"Conflict with -ts_trjaectory_max_cps_ram or TSTrajectorySetMaxCpsRAM. You can set max_cps_ram or max_units_ram, but not both at the same time.");
1886   tjsch->max_units_ram = max_units_ram;
1887   PetscFunctionReturn(0);
1888 }
1889 
1890 static PetscErrorCode TSTrajectorySetMaxUnitsDisk_Memory(TSTrajectory tj,PetscInt max_units_disk)
1891 {
1892   TJScheduler *tjsch = (TJScheduler*)tj->data;
1893 
1894   PetscFunctionBegin;
1895   PetscCheck(tjsch->max_cps_disk,PetscObjectComm((PetscObject)tj),PETSC_ERR_ARG_INCOMP,"Conflict with -ts_trjaectory_max_cps_disk or TSTrajectorySetMaxCpsDisk. You can set max_cps_disk or max_units_disk, but not both at the same time.");
1896   tjsch->max_units_ram = max_units_disk;
1897   PetscFunctionReturn(0);
1898 }
1899 
1900 static PetscErrorCode TSTrajectoryMemorySetType_Memory(TSTrajectory tj,TSTrajectoryMemoryType tj_memory_type)
1901 {
1902   TJScheduler *tjsch = (TJScheduler*)tj->data;
1903 
1904   PetscFunctionBegin;
1905   PetscCheck(!tj->setupcalled,PetscObjectComm((PetscObject)tj),PETSC_ERR_ARG_WRONGSTATE,"Cannot change schedule software after TSTrajectory has been setup or used");
1906   tjsch->tj_memory_type = tj_memory_type;
1907   PetscFunctionReturn(0);
1908 }
1909 
1910 #if defined(PETSC_HAVE_REVOLVE)
1911 PETSC_UNUSED static PetscErrorCode TSTrajectorySetRevolveOnline(TSTrajectory tj,PetscBool use_online)
1912 {
1913   TJScheduler *tjsch = (TJScheduler*)tj->data;
1914 
1915   PetscFunctionBegin;
1916   tjsch->use_online = use_online;
1917   PetscFunctionReturn(0);
1918 }
1919 #endif
1920 
1921 PETSC_UNUSED static PetscErrorCode TSTrajectorySetSaveStack(TSTrajectory tj,PetscBool save_stack)
1922 {
1923   TJScheduler *tjsch = (TJScheduler*)tj->data;
1924 
1925   PetscFunctionBegin;
1926   tjsch->save_stack = save_stack;
1927   PetscFunctionReturn(0);
1928 }
1929 
1930 PETSC_UNUSED static PetscErrorCode TSTrajectorySetUseDRAM(TSTrajectory tj,PetscBool use_dram)
1931 {
1932   TJScheduler *tjsch = (TJScheduler*)tj->data;
1933 
1934   PetscFunctionBegin;
1935   tjsch->stack.use_dram = use_dram;
1936   PetscFunctionReturn(0);
1937 }
1938 
1939 /*@C
1940    TSTrajectoryMemorySetType - sets the software that is used to generate the checkpointing schedule.
1941 
1942    Logically Collective on TSTrajectory
1943 
1944    Input Parameters:
1945 +  tj - the TSTrajectory context
1946 -  tj_memory_type - Revolve or CAMS
1947 
1948    Options Database Key:
1949 .  -ts_trajectory_memory_type <tj_memory_type> - petsc, revolve, cams
1950 
1951    Level: intermediate
1952 
1953    Note:
1954      By default this will use Revolve if it exists
1955 @*/
1956 PetscErrorCode TSTrajectoryMemorySetType(TSTrajectory tj,TSTrajectoryMemoryType tj_memory_type)
1957 {
1958   PetscFunctionBegin;
1959   PetscTryMethod(tj,"TSTrajectoryMemorySetType_C",(TSTrajectory,TSTrajectoryMemoryType),(tj,tj_memory_type));
1960   PetscFunctionReturn(0);
1961 }
1962 
1963 /*@C
1964   TSTrajectorySetMaxCpsRAM - Set maximum number of checkpoints in RAM
1965 
1966   Logically collective
1967 
1968   Input Parameter:
1969 .  tj - tstrajectory context
1970 
1971   Output Parameter:
1972 .  max_cps_ram - maximum number of checkpoints in RAM
1973 
1974   Level: intermediate
1975 
1976 .seealso: `TSTrajectorySetMaxUnitsRAM()`
1977 @*/
1978 PetscErrorCode TSTrajectorySetMaxCpsRAM(TSTrajectory tj,PetscInt max_cps_ram)
1979 {
1980   PetscFunctionBegin;
1981   PetscUseMethod(tj,"TSTrajectorySetMaxCpsRAM_C",(TSTrajectory,PetscInt),(tj,max_cps_ram));
1982   PetscFunctionReturn(0);
1983 }
1984 
1985 /*@C
1986   TSTrajectorySetMaxCpsDisk - Set maximum number of checkpoints on disk
1987 
1988   Logically collective
1989 
1990   Input Parameter:
1991 .  tj - tstrajectory context
1992 
1993   Output Parameter:
1994 .  max_cps_disk - maximum number of checkpoints on disk
1995 
1996   Level: intermediate
1997 
1998 .seealso: `TSTrajectorySetMaxUnitsDisk()`, `TSTrajectorySetMaxUnitsRAM()`
1999 @*/
2000 PetscErrorCode TSTrajectorySetMaxCpsDisk(TSTrajectory tj,PetscInt max_cps_disk)
2001 {
2002   PetscFunctionBegin;
2003   PetscUseMethod(tj,"TSTrajectorySetMaxCpsDisk_C",(TSTrajectory,PetscInt),(tj,max_cps_disk));
2004   PetscFunctionReturn(0);
2005 }
2006 
2007 /*@C
2008   TSTrajectorySetMaxUnitsRAM - Set maximum number of checkpointing units in RAM
2009 
2010   Logically collective
2011 
2012   Input Parameter:
2013 .  tj - tstrajectory context
2014 
2015   Output Parameter:
2016 .  max_units_ram - maximum number of checkpointing units in RAM
2017 
2018   Level: intermediate
2019 
2020 .seealso: `TSTrajectorySetMaxCpsRAM()`
2021 @*/
2022 PetscErrorCode TSTrajectorySetMaxUnitsRAM(TSTrajectory tj,PetscInt max_units_ram)
2023 {
2024   PetscFunctionBegin;
2025   PetscUseMethod(tj,"TSTrajectorySetMaxUnitsRAM_C",(TSTrajectory,PetscInt),(tj,max_units_ram));
2026   PetscFunctionReturn(0);
2027 }
2028 
2029 /*@C
2030   TSTrajectorySetMaxUnitsDisk - Set maximum number of checkpointing units on disk
2031 
2032   Logically collective
2033 
2034   Input Parameter:
2035 .  tj - tstrajectory context
2036 
2037   Output Parameter:
2038 .  max_units_disk - maximum number of checkpointing units on disk
2039 
2040   Level: intermediate
2041 
2042 .seealso: `TSTrajectorySetMaxCpsDisk()`
2043 @*/
2044 PetscErrorCode TSTrajectorySetMaxUnitsDisk(TSTrajectory tj,PetscInt max_units_disk)
2045 {
2046   PetscFunctionBegin;
2047   PetscUseMethod(tj,"TSTrajectorySetMaxUnitsDisk_C",(TSTrajectory,PetscInt),(tj,max_units_disk));
2048   PetscFunctionReturn(0);
2049 }
2050 
2051 static PetscErrorCode TSTrajectorySetFromOptions_Memory(PetscOptionItems *PetscOptionsObject,TSTrajectory tj)
2052 {
2053   TJScheduler    *tjsch = (TJScheduler*)tj->data;
2054   PetscEnum      etmp;
2055   PetscInt       max_cps_ram,max_cps_disk,max_units_ram,max_units_disk;
2056   PetscBool      flg;
2057 
2058   PetscFunctionBegin;
2059   PetscOptionsHeadBegin(PetscOptionsObject,"Memory based TS trajectory options");
2060   {
2061     PetscCall(PetscOptionsInt("-ts_trajectory_max_cps_ram","Maximum number of checkpoints in RAM","TSTrajectorySetMaxCpsRAM",tjsch->max_cps_ram,&max_cps_ram,&flg));
2062     if (flg) {
2063       PetscCall(TSTrajectorySetMaxCpsRAM(tj,max_cps_ram));
2064     }
2065     PetscCall(PetscOptionsInt("-ts_trajectory_max_cps_disk","Maximum number of checkpoints on disk","TSTrajectorySetMaxCpsDisk",tjsch->max_cps_disk,&max_cps_disk,&flg));
2066     if (flg) {
2067       PetscCall(TSTrajectorySetMaxCpsDisk(tj,max_cps_disk));
2068     }
2069     PetscCall(PetscOptionsInt("-ts_trajectory_max_units_ram","Maximum number of checkpointing units in RAM","TSTrajectorySetMaxUnitsRAM",tjsch->max_units_ram,&max_units_ram,&flg));
2070     if (flg) {
2071       PetscCall(TSTrajectorySetMaxUnitsRAM(tj,max_units_ram));
2072     }
2073     PetscCall(PetscOptionsInt("-ts_trajectory_max_units_disk","Maximum number of checkpointing units on disk","TSTrajectorySetMaxUnitsDisk",tjsch->max_units_disk,&max_units_disk,&flg));
2074     if (flg) {
2075       PetscCall(TSTrajectorySetMaxUnitsDisk(tj,max_units_disk));
2076     }
2077     PetscCall(PetscOptionsInt("-ts_trajectory_stride","Stride to save checkpoints to file","TSTrajectorySetStride",tjsch->stride,&tjsch->stride,NULL));
2078 #if defined(PETSC_HAVE_REVOLVE)
2079     PetscCall(PetscOptionsBool("-ts_trajectory_revolve_online","Trick TS trajectory into using online mode of revolve","TSTrajectorySetRevolveOnline",tjsch->use_online,&tjsch->use_online,NULL));
2080 #endif
2081     PetscCall(PetscOptionsBool("-ts_trajectory_save_stack","Save all stack to disk","TSTrajectorySetSaveStack",tjsch->save_stack,&tjsch->save_stack,NULL));
2082     PetscCall(PetscOptionsBool("-ts_trajectory_use_dram","Use DRAM for checkpointing","TSTrajectorySetUseDRAM",tjsch->stack.use_dram,&tjsch->stack.use_dram,NULL));
2083     PetscCall(PetscOptionsEnum("-ts_trajectory_memory_type","Checkpointing scchedule software to use","TSTrajectoryMemorySetType",TSTrajectoryMemoryTypes,(PetscEnum)(int)(tjsch->tj_memory_type),&etmp,&flg));
2084     if (flg) {
2085       PetscCall(TSTrajectoryMemorySetType(tj,(TSTrajectoryMemoryType)etmp));
2086     }
2087   }
2088   PetscOptionsHeadEnd();
2089   PetscFunctionReturn(0);
2090 }
2091 
2092 static PetscErrorCode TSTrajectorySetUp_Memory(TSTrajectory tj,TS ts)
2093 {
2094   TJScheduler    *tjsch = (TJScheduler*)tj->data;
2095   Stack          *stack = &tjsch->stack;
2096 #if defined(PETSC_HAVE_REVOLVE)
2097   RevolveCTX     *rctx,*rctx2;
2098   DiskStack      *diskstack = &tjsch->diskstack;
2099   PetscInt       diskblocks;
2100 #endif
2101   PetscInt       numY,total_steps;
2102   PetscBool      fixedtimestep;
2103 
2104   PetscFunctionBegin;
2105   if (ts->adapt) {
2106     PetscCall(PetscObjectTypeCompare((PetscObject)ts->adapt,TSADAPTNONE,&fixedtimestep));
2107   } else {
2108     fixedtimestep = PETSC_TRUE;
2109   }
2110   total_steps = (PetscInt)(PetscCeilReal((ts->max_time-ts->ptime)/ts->time_step));
2111   total_steps = total_steps < 0 ? PETSC_MAX_INT : total_steps;
2112   if (fixedtimestep) tjsch->total_steps = PetscMin(ts->max_steps,total_steps);
2113 
2114   tjsch->stack.solution_only = tj->solution_only;
2115   PetscCall(TSGetStages(ts,&numY,PETSC_IGNORE));
2116   if (stack->solution_only) {
2117     if (tjsch->max_units_ram) tjsch->max_cps_ram = tjsch->max_units_ram;
2118     else tjsch->max_units_ram = tjsch->max_cps_ram;
2119     if (tjsch->max_units_disk) tjsch->max_cps_disk = tjsch->max_units_disk;
2120   } else {
2121     if (tjsch->max_units_ram) tjsch->max_cps_ram = (ts->stifflyaccurate) ? tjsch->max_units_ram/numY : tjsch->max_units_ram/(numY+1);
2122     else tjsch->max_units_ram = (ts->stifflyaccurate) ? numY*tjsch->max_cps_ram : (numY+1)*tjsch->max_cps_ram;
2123     if (tjsch->max_units_disk) tjsch->max_cps_disk = (ts->stifflyaccurate) ? tjsch->max_units_disk/numY : tjsch->max_units_disk/(numY+1);
2124     else tjsch->max_units_disk = (ts->stifflyaccurate) ? numY*tjsch->max_cps_disk : (numY+1)*tjsch->max_cps_disk;
2125   }
2126   if (tjsch->max_cps_ram > 0) stack->stacksize = tjsch->max_units_ram; /* maximum stack size. Could be overallocated. */
2127 
2128   /* Determine the scheduler type */
2129   if (tjsch->stride > 1) { /* two level mode */
2130     PetscCheck(!tjsch->save_stack || tjsch->max_cps_disk <= 1 || tjsch->max_cps_disk > tjsch->max_cps_ram,PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_INCOMP,"The specified disk capacity is not enough to store a full stack of RAM checkpoints. You might want to change the disk capacity or use single level checkpointing instead.");
2131     if (tjsch->max_cps_disk <= 1 && tjsch->max_cps_ram > 1 && tjsch->max_cps_ram <= tjsch->stride-1) tjsch->stype = TWO_LEVEL_REVOLVE; /* use revolve_offline for each stride */
2132     if (tjsch->max_cps_disk > 1 && tjsch->max_cps_ram > 1 && tjsch->max_cps_ram <= tjsch->stride-1) tjsch->stype = TWO_LEVEL_TWO_REVOLVE;  /* use revolve_offline for each stride */
2133     if (tjsch->max_cps_disk <= 1 && (tjsch->max_cps_ram >= tjsch->stride || tjsch->max_cps_ram == -1)) tjsch->stype = TWO_LEVEL_NOREVOLVE; /* can also be handled by TWO_LEVEL_REVOLVE */
2134   } else { /* single level mode */
2135     if (fixedtimestep) {
2136       if (tjsch->max_cps_ram >= tjsch->total_steps-1 || tjsch->max_cps_ram == -1)
2137         tjsch->stype = NONE; /* checkpoint all */
2138       else { /* choose the schedule software for offline checkpointing */
2139         switch (tjsch->tj_memory_type) {
2140           case TJ_PETSC:
2141             tjsch->stype = NONE;
2142             break;
2143           case TJ_CAMS:
2144             tjsch->stype = CAMS_OFFLINE;
2145             break;
2146           case TJ_REVOLVE:
2147             tjsch->stype = (tjsch->max_cps_disk>1) ? REVOLVE_MULTISTAGE : REVOLVE_OFFLINE;
2148             break;
2149           default:
2150             break;
2151         }
2152       }
2153     } else tjsch->stype = NONE; /* checkpoint all for adaptive time step */
2154 #if defined(PETSC_HAVE_REVOLVE)
2155     if (tjsch->use_online) tjsch->stype = REVOLVE_ONLINE; /* trick into online (for testing purpose only) */
2156 #endif
2157     PetscCheck(tjsch->stype == NONE || tjsch->max_cps_ram >= 1 || tjsch->max_cps_disk >= 1,PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_INCOMP,"The specified storage capacity is insufficient for one checkpoint, which is the minimum");
2158   }
2159   if (tjsch->stype >= CAMS_OFFLINE) {
2160 #ifndef PETSC_HAVE_CAMS
2161     SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"CAMS is needed when there is not enough memory to checkpoint all time steps according to the user's settings, please reconfigure with the additional option --download-cams.");
2162 #else
2163     CAMSCTX  *actx;
2164     PetscInt ns = 0;
2165     if (stack->solution_only) {
2166       offline_ca_create(tjsch->total_steps,tjsch->max_cps_ram);
2167     } else {
2168       PetscCall(TSGetStages(ts,&ns,PETSC_IGNORE));
2169       offline_cams_create(tjsch->total_steps,tjsch->max_units_ram,ns,ts->stifflyaccurate);
2170     }
2171     PetscCall(PetscNew(&actx));
2172     actx->lastcheckpointstep    = -1; /* -1 can trigger the initialization of CAMS */
2173     actx->lastcheckpointtype    = -1; /* -1 can trigger the initialization of CAMS */
2174     actx->endstep               = tjsch->total_steps;
2175     actx->num_units_avail       = tjsch->max_units_ram;
2176     actx->num_stages            = ns;
2177     tjsch->actx                 = actx;
2178 #endif
2179   } else if (tjsch->stype > TWO_LEVEL_NOREVOLVE) {
2180 #ifndef PETSC_HAVE_REVOLVE
2181     SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"revolve is needed when there is not enough memory to checkpoint all time steps according to the user's settings, please reconfigure with the additional option --download-revolve.");
2182 #else
2183     PetscRevolveInt rfine,rsnaps,rsnaps2;
2184 
2185     switch (tjsch->stype) {
2186       case TWO_LEVEL_REVOLVE:
2187         PetscCall(PetscRevolveIntCast(tjsch->stride,&rfine));
2188         PetscCall(PetscRevolveIntCast(tjsch->max_cps_ram,&rsnaps));
2189         revolve_create_offline(rfine,rsnaps);
2190         break;
2191       case TWO_LEVEL_TWO_REVOLVE:
2192         diskblocks = tjsch->save_stack ? tjsch->max_cps_disk/(tjsch->max_cps_ram+1) : tjsch->max_cps_disk; /* The block size depends on whether the stack is saved. */
2193         diskstack->stacksize = diskblocks;
2194         PetscCall(PetscRevolveIntCast(tjsch->stride,&rfine));
2195         PetscCall(PetscRevolveIntCast(tjsch->max_cps_ram,&rsnaps));
2196         revolve_create_offline(rfine,rsnaps);
2197         PetscCall(PetscRevolveIntCast((tjsch->total_steps+tjsch->stride-1)/tjsch->stride,&rfine));
2198         PetscCall(PetscRevolveIntCast(diskblocks,&rsnaps));
2199         revolve2_create_offline(rfine,rsnaps);
2200         PetscCall(PetscNew(&rctx2));
2201         rctx2->snaps_in       = rsnaps;
2202         rctx2->reverseonestep = PETSC_FALSE;
2203         rctx2->check          = 0;
2204         rctx2->oldcapo        = 0;
2205         rctx2->capo           = 0;
2206         rctx2->info           = 2;
2207         rctx2->fine           = rfine;
2208         tjsch->rctx2          = rctx2;
2209         diskstack->top        = -1;
2210         PetscCall(PetscMalloc1(diskstack->stacksize,&diskstack->container));
2211         break;
2212       case REVOLVE_OFFLINE:
2213         PetscCall(PetscRevolveIntCast(tjsch->total_steps,&rfine));
2214         PetscCall(PetscRevolveIntCast(tjsch->max_cps_ram,&rsnaps));
2215         revolve_create_offline(rfine,rsnaps);
2216         break;
2217       case REVOLVE_ONLINE:
2218         stack->stacksize = tjsch->max_cps_ram;
2219         PetscCall(PetscRevolveIntCast(tjsch->max_cps_ram,&rsnaps));
2220         revolve_create_online(rsnaps);
2221         break;
2222       case REVOLVE_MULTISTAGE:
2223         PetscCall(PetscRevolveIntCast(tjsch->total_steps,&rfine));
2224         PetscCall(PetscRevolveIntCast(tjsch->max_cps_ram,&rsnaps));
2225         PetscCall(PetscRevolveIntCast(tjsch->max_cps_ram+tjsch->max_cps_disk,&rsnaps2));
2226         revolve_create_multistage(rfine,rsnaps2,rsnaps);
2227         break;
2228       default:
2229         break;
2230     }
2231     PetscCall(PetscNew(&rctx));
2232     PetscCall(PetscRevolveIntCast(tjsch->max_cps_ram,&rsnaps));
2233     rctx->snaps_in       = rsnaps; /* for theta methods snaps_in=2*max_cps_ram */
2234     rctx->reverseonestep = PETSC_FALSE;
2235     rctx->check          = 0;
2236     rctx->oldcapo        = 0;
2237     rctx->capo           = 0;
2238     rctx->info           = 2;
2239     if (tjsch->stride > 1) {
2240       PetscCall(PetscRevolveIntCast(tjsch->stride,&rfine));
2241     } else {
2242       PetscCall(PetscRevolveIntCast(tjsch->total_steps,&rfine));
2243     }
2244     rctx->fine           = rfine;
2245     tjsch->rctx          = rctx;
2246     if (tjsch->stype == REVOLVE_ONLINE) rctx->fine = -1;
2247 #endif
2248   } else {
2249     if (tjsch->stype == TWO_LEVEL_NOREVOLVE) stack->stacksize = tjsch->stride-1; /* need tjsch->stride-1 at most */
2250     if (tjsch->stype == NONE) {
2251       if (fixedtimestep) stack->stacksize = stack->solution_only ? tjsch->total_steps : tjsch->total_steps-1;
2252       else { /* adaptive time step */
2253         /* if max_cps_ram is not specified, use maximal allowed number of steps for stack size */
2254         if (tjsch->max_cps_ram == -1) stack->stacksize = ts->max_steps < PETSC_MAX_INT ? ts->max_steps : 10000;
2255         tjsch->total_steps = stack->solution_only ? stack->stacksize : stack->stacksize+1; /* will be updated as time integration advances */
2256       }
2257     }
2258   }
2259 
2260   if ((tjsch->stype >= TWO_LEVEL_NOREVOLVE && tjsch->stype < REVOLVE_OFFLINE) || tjsch->stype == REVOLVE_MULTISTAGE) { /* these types need to use disk */
2261     PetscCall(TSTrajectorySetUp_Basic(tj,ts));
2262   }
2263 
2264   stack->stacksize = PetscMax(stack->stacksize,1);
2265   tjsch->recompute = PETSC_FALSE;
2266   PetscCall(StackInit(stack,stack->stacksize,numY));
2267   PetscFunctionReturn(0);
2268 }
2269 
2270 static PetscErrorCode TSTrajectoryReset_Memory(TSTrajectory tj)
2271 {
2272 #if defined (PETSC_HAVE_REVOLVE) || defined (PETSC_HAVE_CAMS)
2273   TJScheduler *tjsch = (TJScheduler*)tj->data;
2274 #endif
2275 
2276   PetscFunctionBegin;
2277 #if defined(PETSC_HAVE_REVOLVE)
2278   if (tjsch->stype > TWO_LEVEL_NOREVOLVE) {
2279     revolve_reset();
2280     if (tjsch->stype == TWO_LEVEL_TWO_REVOLVE) {
2281       revolve2_reset();
2282       PetscCall(PetscFree(tjsch->diskstack.container));
2283     }
2284   }
2285   if (tjsch->stype > TWO_LEVEL_NOREVOLVE) {
2286     PetscCall(PetscFree(tjsch->rctx));
2287     PetscCall(PetscFree(tjsch->rctx2));
2288   }
2289 #endif
2290 #if defined(PETSC_HAVE_CAMS)
2291   if (tjsch->stype == CAMS_OFFLINE) {
2292     if (tjsch->stack.solution_only) offline_ca_destroy();
2293     else offline_ca_destroy();
2294     PetscCall(PetscFree(tjsch->actx));
2295   }
2296 #endif
2297   PetscFunctionReturn(0);
2298 }
2299 
2300 static PetscErrorCode TSTrajectoryDestroy_Memory(TSTrajectory tj)
2301 {
2302   TJScheduler    *tjsch = (TJScheduler*)tj->data;
2303 
2304   PetscFunctionBegin;
2305   PetscCall(StackDestroy(&tjsch->stack));
2306   PetscCall(PetscViewerDestroy(&tjsch->viewer));
2307   PetscCall(PetscObjectComposeFunction((PetscObject)tj,"TSTrajectorySetMaxCpsRAM_C",NULL));
2308   PetscCall(PetscObjectComposeFunction((PetscObject)tj,"TSTrajectorySetMaxCpsDisk_C",NULL));
2309   PetscCall(PetscObjectComposeFunction((PetscObject)tj,"TSTrajectorySetMaxUnitsRAM_C",NULL));
2310   PetscCall(PetscObjectComposeFunction((PetscObject)tj,"TSTrajectorySetMaxUnitsDisk_C",NULL));
2311   PetscCall(PetscObjectComposeFunction((PetscObject)tj,"TSTrajectoryMemorySetType_C",NULL));
2312   PetscCall(PetscFree(tjsch));
2313   PetscFunctionReturn(0);
2314 }
2315 
2316 /*MC
2317       TSTRAJECTORYMEMORY - Stores each solution of the ODE/ADE in memory
2318 
2319   Level: intermediate
2320 
2321 .seealso: `TSTrajectoryCreate()`, `TS`, `TSTrajectorySetType()`
2322 
2323 M*/
2324 PETSC_EXTERN PetscErrorCode TSTrajectoryCreate_Memory(TSTrajectory tj,TS ts)
2325 {
2326   TJScheduler    *tjsch;
2327 
2328   PetscFunctionBegin;
2329   tj->ops->set            = TSTrajectorySet_Memory;
2330   tj->ops->get            = TSTrajectoryGet_Memory;
2331   tj->ops->setup          = TSTrajectorySetUp_Memory;
2332   tj->ops->setfromoptions = TSTrajectorySetFromOptions_Memory;
2333   tj->ops->reset          = TSTrajectoryReset_Memory;
2334   tj->ops->destroy        = TSTrajectoryDestroy_Memory;
2335 
2336   PetscCall(PetscNew(&tjsch));
2337   tjsch->stype        = NONE;
2338   tjsch->max_cps_ram  = -1; /* -1 indicates that it is not set */
2339   tjsch->max_cps_disk = -1; /* -1 indicates that it is not set */
2340   tjsch->stride       = 0; /* if not zero, two-level checkpointing will be used */
2341 #if defined(PETSC_HAVE_REVOLVE)
2342   tjsch->use_online   = PETSC_FALSE;
2343 #endif
2344   tjsch->save_stack   = PETSC_TRUE;
2345 
2346   tjsch->stack.solution_only = tj->solution_only;
2347   PetscCall(PetscViewerCreate(PetscObjectComm((PetscObject)tj),&tjsch->viewer));
2348   PetscCall(PetscViewerSetType(tjsch->viewer,PETSCVIEWERBINARY));
2349   PetscCall(PetscViewerPushFormat(tjsch->viewer,PETSC_VIEWER_NATIVE));
2350   PetscCall(PetscViewerFileSetMode(tjsch->viewer,FILE_MODE_WRITE));
2351 
2352   PetscCall(PetscObjectComposeFunction((PetscObject)tj,"TSTrajectorySetMaxCpsRAM_C",TSTrajectorySetMaxCpsRAM_Memory));
2353   PetscCall(PetscObjectComposeFunction((PetscObject)tj,"TSTrajectorySetMaxCpsDisk_C",TSTrajectorySetMaxCpsDisk_Memory));
2354   PetscCall(PetscObjectComposeFunction((PetscObject)tj,"TSTrajectorySetMaxUnitsRAM_C",TSTrajectorySetMaxUnitsRAM_Memory));
2355   PetscCall(PetscObjectComposeFunction((PetscObject)tj,"TSTrajectorySetMaxUnitsDisk_C",TSTrajectorySetMaxUnitsDisk_Memory));
2356   PetscCall(PetscObjectComposeFunction((PetscObject)tj,"TSTrajectoryMemorySetType_C",TSTrajectoryMemorySetType_Memory));
2357   tj->data = tjsch;
2358   PetscFunctionReturn(0);
2359 }
2360