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