xref: /petsc/src/ts/trajectory/impls/memory/trajmemory.c (revision 13850c04a7377a91d5b1bff18e84c75ea74ba331)
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 #endif
6 #if defined(PETSC_HAVE_MEMKIND)
7 extern PetscErrorCode PetscMallocSetDRAM(void);
8 extern PetscErrorCode PetscMallocResetDRAM(void);
9 #endif
10 
11 PetscLogEvent TSTrajectory_DiskWrite, TSTrajectory_DiskRead;
12 
13 typedef enum {NONE,TWO_LEVEL_NOREVOLVE,TWO_LEVEL_REVOLVE,TWO_LEVEL_TWO_REVOLVE,REVOLVE_OFFLINE,REVOLVE_ONLINE,REVOLVE_MULTISTAGE} SchedulerType;
14 
15 typedef struct _StackElement {
16   PetscInt  stepnum;
17   Vec       X;
18   Vec       *Y;
19   PetscReal time;
20   PetscReal timeprev; /* for no solution_only mode */
21   PetscReal timenext; /* for solution_only mode */
22 } *StackElement;
23 
24 #if defined(PETSC_HAVE_REVOLVE)
25 typedef struct _RevolveCTX {
26   PetscBool reverseonestep;
27   PetscInt  where;
28   PetscInt  snaps_in;
29   PetscInt  stepsleft;
30   PetscInt  check;
31   PetscInt  oldcapo;
32   PetscInt  capo;
33   PetscInt  fine;
34   PetscInt  info;
35 } RevolveCTX;
36 #endif
37 
38 typedef struct _Stack {
39   PetscInt      stacksize;
40   PetscInt      top;
41   StackElement  *container;
42   PetscInt      numY;
43   PetscBool     solution_only;
44 } Stack;
45 
46 typedef struct _DiskStack {
47   PetscInt  stacksize;
48   PetscInt  top;
49   PetscInt  *container;
50 } DiskStack;
51 
52 typedef struct _TJScheduler {
53   SchedulerType stype;
54 #if defined(PETSC_HAVE_REVOLVE)
55   RevolveCTX    *rctx,*rctx2;
56   PetscBool     use_online;
57   PetscInt      store_stride;
58 #endif
59   PetscBool     recompute;
60   PetscBool     skip_trajectory;
61   PetscBool     save_stack;
62   PetscInt      max_cps_ram;  /* maximum checkpoints in RAM */
63   PetscInt      max_cps_disk; /* maximum checkpoints on disk */
64   PetscInt      stride;
65   PetscInt      total_steps;  /* total number of steps */
66   Stack         stack;
67   DiskStack     diskstack;
68 } TJScheduler;
69 
70 static PetscErrorCode TurnForwardWithStepsize(TS ts,PetscReal nextstepsize)
71 {
72     PetscReal      stepsize;
73     PetscErrorCode ierr;
74 
75     PetscFunctionBegin;
76     /* reverse the direction */
77     ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr);
78     stepsize = nextstepsize;
79     ierr = TSSetTimeStep(ts,stepsize);CHKERRQ(ierr);
80     PetscFunctionReturn(0);
81 }
82 
83 static PetscErrorCode TurnForward(TS ts)
84 {
85   PetscReal      stepsize;
86   PetscErrorCode ierr;
87 
88   PetscFunctionBegin;
89   /* reverse the direction */
90   ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr);
91   ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr);
92   PetscFunctionReturn(0);
93 }
94 
95 static PetscErrorCode TurnBackward(TS ts)
96 {
97   PetscReal      stepsize;
98   PetscErrorCode ierr;
99 
100   PetscFunctionBegin;
101   /* reverse the direction */
102   stepsize = ts->ptime_prev-ts->ptime;
103   ierr = TSSetTimeStep(ts,stepsize);CHKERRQ(ierr);
104   PetscFunctionReturn(0);
105 }
106 
107 static PetscErrorCode ElementCreate(TS ts,Stack *stack,StackElement *e)
108 {
109   Vec            X;
110   Vec            *Y;
111   PetscErrorCode ierr;
112 
113   PetscFunctionBegin;
114 #if defined(PETSC_HAVE_MEMKIND)
115   ierr = PetscMallocSetDRAM();CHKERRQ(ierr);
116 #endif
117   ierr = PetscCalloc1(1,e);CHKERRQ(ierr);
118   ierr = TSGetSolution(ts,&X);CHKERRQ(ierr);
119   ierr = VecDuplicate(X,&(*e)->X);CHKERRQ(ierr);
120   if (stack->numY > 0 && !stack->solution_only) {
121     ierr = TSGetStages(ts,&stack->numY,&Y);CHKERRQ(ierr);
122     ierr = VecDuplicateVecs(Y[0],stack->numY,&(*e)->Y);CHKERRQ(ierr);
123   }
124 #if defined(PETSC_HAVE_MEMKIND)
125   ierr = PetscMallocResetDRAM();CHKERRQ(ierr);
126 #endif
127   PetscFunctionReturn(0);
128 }
129 
130 static PetscErrorCode ElementSet(TS ts,Stack *stack,StackElement *e,PetscInt stepnum,PetscReal time,Vec X)
131 {
132   Vec            *Y;
133   PetscInt       i;
134   PetscReal      timeprev;
135   PetscErrorCode ierr;
136 
137   PetscFunctionBegin;
138   ierr = VecCopy(X,(*e)->X);CHKERRQ(ierr);
139   if (stack->numY > 0 && !stack->solution_only) {
140     ierr = TSGetStages(ts,&stack->numY,&Y);CHKERRQ(ierr);
141     for (i=0;i<stack->numY;i++) {
142       ierr = VecCopy(Y[i],(*e)->Y[i]);CHKERRQ(ierr);
143     }
144   }
145   (*e)->stepnum = stepnum;
146   (*e)->time    = time;
147   /* for consistency */
148   if (stepnum == 0) {
149     (*e)->timeprev = (*e)->time - ts->time_step;
150   } else {
151     ierr = TSGetPrevTime(ts,&timeprev);CHKERRQ(ierr);
152     (*e)->timeprev = timeprev;
153   }
154   PetscFunctionReturn(0);
155 }
156 
157 static PetscErrorCode ElementDestroy(Stack *stack,StackElement e)
158 {
159   PetscErrorCode ierr;
160 
161   PetscFunctionBegin;
162 #if defined(PETSC_HAVE_MEMKIND)
163   ierr = PetscMallocSetDRAM();CHKERRQ(ierr);
164 #endif
165   ierr = VecDestroy(&e->X);CHKERRQ(ierr);
166   if (stack->numY > 0 && !stack->solution_only) {
167     ierr = VecDestroyVecs(stack->numY,&e->Y);CHKERRQ(ierr);
168   }
169   ierr = PetscFree(e);CHKERRQ(ierr);
170 #if defined(PETSC_HAVE_MEMKIND)
171   ierr = PetscMallocResetDRAM();CHKERRQ(ierr);
172 #endif
173   PetscFunctionReturn(0);
174 }
175 
176 static PetscErrorCode StackResize(Stack *stack,PetscInt newsize)
177 {
178   StackElement   *newcontainer;
179   PetscInt       i;
180   PetscErrorCode ierr;
181 
182   PetscFunctionBegin;
183   ierr = PetscMalloc1(newsize*sizeof(StackElement),&newcontainer);CHKERRQ(ierr);
184   for (i=0;i<stack->stacksize;i++) {
185     newcontainer[i] = stack->container[i];
186   }
187   ierr = PetscFree(stack->container);CHKERRQ(ierr);
188   stack->container = newcontainer;
189   stack->stacksize = newsize;
190   PetscFunctionReturn(0);
191 }
192 
193 static PetscErrorCode StackPush(Stack *stack,StackElement e)
194 {
195   PetscFunctionBegin;
196   if (stack->top+1 >= stack->stacksize) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MEMC,"Maximum stack size (%D) exceeded",stack->stacksize);
197   stack->container[++stack->top] = e;
198   PetscFunctionReturn(0);
199 }
200 
201 static PetscErrorCode StackPop(Stack *stack,StackElement *e)
202 {
203   PetscFunctionBegin;
204   if (stack->top == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEMC,"Empty stack");
205   *e = stack->container[stack->top--];
206   PetscFunctionReturn(0);
207 }
208 
209 static PetscErrorCode StackTop(Stack *stack,StackElement *e)
210 {
211   PetscFunctionBegin;
212   *e = stack->container[stack->top];
213   PetscFunctionReturn(0);
214 }
215 
216 static PetscErrorCode StackCreate(Stack *stack,PetscInt size,PetscInt ny)
217 {
218   PetscErrorCode ierr;
219 
220   PetscFunctionBegin;
221   stack->top  = -1;
222   stack->numY = ny;
223 
224   ierr = PetscMalloc1(size*sizeof(StackElement),&stack->container);CHKERRQ(ierr);
225   PetscFunctionReturn(0);
226 }
227 
228 static PetscErrorCode StackDestroy(Stack *stack)
229 {
230   PetscInt       i,n;
231   StackElement   e;
232   PetscErrorCode ierr;
233 
234   PetscFunctionBegin;
235   if (stack->top > -1) {
236     n = stack->top+1; /* number of elements in the stack */
237     for (i=0;i<n;i++) {
238       ierr = StackPop(stack,&e);CHKERRQ(ierr);
239       ierr = ElementDestroy(stack,e);CHKERRQ(ierr);
240     }
241   }
242   ierr = PetscFree(stack->container);CHKERRQ(ierr);
243   PetscFunctionReturn(0);
244 }
245 
246 #if defined(PETSC_HAVE_REVOLVE)
247 static PetscErrorCode StackFind(Stack *stack,StackElement *e,PetscInt index)
248 {
249   PetscFunctionBegin;
250   *e = stack->container[index];
251   PetscFunctionReturn(0);
252 }
253 #endif
254 
255 static PetscErrorCode OutputBIN(const char *filename,PetscViewer *viewer)
256 {
257   PetscErrorCode ierr;
258 
259   PetscFunctionBegin;
260   ierr = PetscViewerCreate(PETSC_COMM_WORLD,viewer);CHKERRQ(ierr);
261   ierr = PetscViewerSetType(*viewer,PETSCVIEWERBINARY);CHKERRQ(ierr);
262   ierr = PetscViewerFileSetMode(*viewer,FILE_MODE_WRITE);CHKERRQ(ierr);
263   ierr = PetscViewerFileSetName(*viewer,filename);CHKERRQ(ierr);
264   PetscFunctionReturn(0);
265 }
266 
267 static PetscErrorCode WriteToDisk(PetscInt stepnum,PetscReal time,PetscReal timeprev,Vec X,Vec *Y,PetscInt numY,PetscBool solution_only,PetscViewer viewer)
268 {
269   PetscInt       i;
270   PetscErrorCode ierr;
271 
272   PetscFunctionBegin;
273   ierr = PetscViewerBinaryWrite(viewer,&stepnum,1,PETSC_INT,PETSC_FALSE);CHKERRQ(ierr);
274   ierr = VecView(X,viewer);CHKERRQ(ierr);
275   for (i=0;!solution_only && i<numY;i++) {
276     ierr = VecView(Y[i],viewer);CHKERRQ(ierr);
277   }
278   ierr = PetscViewerBinaryWrite(viewer,&time,1,PETSC_REAL,PETSC_FALSE);CHKERRQ(ierr);
279   ierr = PetscViewerBinaryWrite(viewer,&timeprev,1,PETSC_REAL,PETSC_FALSE);CHKERRQ(ierr);
280   PetscFunctionReturn(0);
281 }
282 
283 static PetscErrorCode ReadFromDisk(PetscInt *stepnum,PetscReal *time,PetscReal *timeprev,Vec X,Vec *Y,PetscInt numY,PetscBool solution_only,PetscViewer viewer)
284 {
285   PetscInt       i;
286   PetscErrorCode ierr;
287 
288   PetscFunctionBegin;
289   ierr = PetscViewerBinaryRead(viewer,stepnum,1,NULL,PETSC_INT);CHKERRQ(ierr);
290   ierr = VecLoad(X,viewer);CHKERRQ(ierr);
291   for (i=0;!solution_only && i<numY;i++) {
292     ierr = VecLoad(Y[i],viewer);CHKERRQ(ierr);
293   }
294   ierr = PetscViewerBinaryRead(viewer,time,1,NULL,PETSC_REAL);CHKERRQ(ierr);
295   ierr = PetscViewerBinaryRead(viewer,timeprev,1,NULL,PETSC_REAL);CHKERRQ(ierr);
296   PetscFunctionReturn(0);
297 }
298 
299 static PetscErrorCode StackDumpAll(TSTrajectory tj,TS ts,Stack *stack,PetscInt id)
300 {
301   Vec            *Y;
302   PetscInt       i;
303   StackElement   e = NULL;
304   PetscViewer    viewer;
305   char           filename[PETSC_MAX_PATH_LEN];
306   PetscErrorCode ierr;
307 
308   PetscFunctionBegin;
309   if (tj->monitor) {
310     ierr = PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr);
311     ierr = PetscViewerASCIIPrintf(tj->monitor,"Dump stack to file\n");CHKERRQ(ierr);
312     ierr = PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr);
313   }
314   if (id == 1) {
315     PetscMPIInt rank;
316     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)ts),&rank);CHKERRQ(ierr);
317     if (!rank) {
318       ierr = PetscRMTree("SA-data");CHKERRQ(ierr);
319       ierr = PetscMkdir("SA-data");CHKERRQ(ierr);
320     }
321   }
322   ierr = PetscSNPrintf(filename,sizeof(filename),"SA-data/SA-STACK%06d.bin",id);CHKERRQ(ierr);
323   ierr = OutputBIN(filename,&viewer);CHKERRQ(ierr);
324   for (i=0;i<stack->stacksize;i++) {
325     e = stack->container[i];
326     ierr = PetscLogEventBegin(TSTrajectory_DiskWrite,ts,0,0,0);CHKERRQ(ierr);
327     ierr = WriteToDisk(e->stepnum,e->time,e->timeprev,e->X,e->Y,stack->numY,stack->solution_only,viewer);CHKERRQ(ierr);
328     ierr = PetscLogEventEnd(TSTrajectory_DiskWrite,ts,0,0,0);CHKERRQ(ierr);
329     ts->trajectory->diskwrites++;
330   }
331   /* 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 */
332   ierr = TSGetStages(ts,&stack->numY,&Y);CHKERRQ(ierr);
333   ierr = PetscLogEventBegin(TSTrajectory_DiskWrite,ts,0,0,0);CHKERRQ(ierr);
334   ierr = WriteToDisk(ts->total_steps,ts->ptime,ts->ptime_prev,ts->vec_sol,Y,stack->numY,stack->solution_only,viewer);CHKERRQ(ierr);
335   ierr = PetscLogEventEnd(TSTrajectory_DiskWrite,ts,0,0,0);CHKERRQ(ierr);
336   ts->trajectory->diskwrites++;
337   for (i=0;i<stack->stacksize;i++) {
338     ierr = StackPop(stack,&e);CHKERRQ(ierr);
339     ierr = ElementDestroy(stack,e);CHKERRQ(ierr);
340   }
341   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
342   PetscFunctionReturn(0);
343 }
344 
345 static PetscErrorCode StackLoadAll(TSTrajectory tj,TS ts,Stack *stack,PetscInt id)
346 {
347   Vec            *Y;
348   PetscInt       i;
349   StackElement   e;
350   PetscViewer    viewer;
351   char           filename[PETSC_MAX_PATH_LEN];
352   PetscErrorCode ierr;
353 
354   PetscFunctionBegin;
355   if (tj->monitor) {
356     ierr = PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr);
357     ierr = PetscViewerASCIIPrintf(tj->monitor,"Load stack from file\n");CHKERRQ(ierr);
358     ierr = PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr);
359   }
360   ierr = PetscSNPrintf(filename,sizeof filename,"SA-data/SA-STACK%06d.bin",id);CHKERRQ(ierr);
361   ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);CHKERRQ(ierr);
362   for (i=0;i<stack->stacksize;i++) {
363     ierr = ElementCreate(ts,stack,&e);CHKERRQ(ierr);
364     ierr = StackPush(stack,e);CHKERRQ(ierr);
365     ierr = PetscLogEventBegin(TSTrajectory_DiskRead,ts,0,0,0);CHKERRQ(ierr);
366     ierr = ReadFromDisk(&e->stepnum,&e->time,&e->timeprev,e->X,e->Y,stack->numY,stack->solution_only,viewer);CHKERRQ(ierr);
367     ierr = PetscLogEventEnd(TSTrajectory_DiskRead,ts,0,0,0);CHKERRQ(ierr);
368     ts->trajectory->diskreads++;
369   }
370   /* load the last step into TS */
371   ierr = TSGetStages(ts,&stack->numY,&Y);CHKERRQ(ierr);
372   ierr = PetscLogEventBegin(TSTrajectory_DiskRead,ts,0,0,0);CHKERRQ(ierr);
373   ierr = ReadFromDisk(&ts->total_steps,&ts->ptime,&ts->ptime_prev,ts->vec_sol,Y,stack->numY,stack->solution_only,viewer);CHKERRQ(ierr);
374   ierr = PetscLogEventEnd(TSTrajectory_DiskRead,ts,0,0,0);CHKERRQ(ierr);
375   ts->trajectory->diskreads++;
376   ierr = TurnBackward(ts);CHKERRQ(ierr);
377   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
378   PetscFunctionReturn(0);
379 }
380 
381 #if defined(PETSC_HAVE_REVOLVE)
382 static PetscErrorCode StackLoadLast(TSTrajectory tj,TS ts,Stack *stack,PetscInt id)
383 {
384   Vec            *Y;
385   PetscInt       size;
386   PetscViewer    viewer;
387   char           filename[PETSC_MAX_PATH_LEN];
388 #if defined(PETSC_HAVE_MPIIO)
389   PetscBool      usempiio;
390 #endif
391   int            fd;
392   off_t          off,offset;
393   PetscErrorCode ierr;
394 
395   PetscFunctionBegin;
396   if (tj->monitor) {
397     ierr = PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr);
398     ierr = PetscViewerASCIIPrintf(tj->monitor,"Load last stack element from file\n");CHKERRQ(ierr);
399     ierr = PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr);
400   }
401   ierr = TSGetStages(ts,&stack->numY,&Y);CHKERRQ(ierr);
402   ierr = VecGetSize(Y[0],&size);CHKERRQ(ierr);
403   /* VecView writes to file two extra int's for class id and number of rows */
404   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;
405 
406   ierr = PetscSNPrintf(filename,sizeof filename,"SA-data/SA-STACK%06d.bin",id);CHKERRQ(ierr);
407   ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);CHKERRQ(ierr);
408 #if defined(PETSC_HAVE_MPIIO)
409   ierr = PetscViewerBinaryGetUseMPIIO(viewer,&usempiio);
410   if (usempiio) {
411     ierr = PetscViewerBinaryGetMPIIODescriptor(viewer,(MPI_File*)&fd);CHKERRQ(ierr);
412     ierr = PetscBinarySynchronizedSeek(PETSC_COMM_WORLD,fd,off,PETSC_BINARY_SEEK_END,&offset);CHKERRQ(ierr);
413   } else {
414 #endif
415     ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
416     ierr = PetscBinarySeek(fd,off,PETSC_BINARY_SEEK_END,&offset);CHKERRQ(ierr);
417 #if defined(PETSC_HAVE_MPIIO)
418   }
419 #endif
420   /* load the last step into TS */
421   ierr = PetscLogEventBegin(TSTrajectory_DiskRead,ts,0,0,0);CHKERRQ(ierr);
422   ierr = ReadFromDisk(&ts->total_steps,&ts->ptime,&ts->ptime_prev,ts->vec_sol,Y,stack->numY,stack->solution_only,viewer);CHKERRQ(ierr);
423   ierr = PetscLogEventEnd(TSTrajectory_DiskRead,ts,0,0,0);CHKERRQ(ierr);
424   ts->trajectory->diskreads++;
425   ierr = TurnBackward(ts);CHKERRQ(ierr);
426   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
427   PetscFunctionReturn(0);
428 }
429 #endif
430 
431 static PetscErrorCode DumpSingle(TSTrajectory tj,TS ts,Stack *stack,PetscInt id)
432 {
433   Vec            *Y;
434   PetscInt       stepnum;
435   PetscViewer    viewer;
436   char           filename[PETSC_MAX_PATH_LEN];
437   PetscErrorCode ierr;
438 
439   PetscFunctionBegin;
440   if (tj->monitor) {
441     ierr = PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr);
442     ierr = PetscViewerASCIIPrintf(tj->monitor,"Load a single point from file\n");CHKERRQ(ierr);
443     ierr = PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr);
444   }
445   ierr = TSGetTotalSteps(ts,&stepnum);CHKERRQ(ierr);
446   if (id == 1) {
447     PetscMPIInt rank;
448     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)ts),&rank);CHKERRQ(ierr);
449     if (!rank) {
450       ierr = PetscRMTree("SA-data");CHKERRQ(ierr);
451       ierr = PetscMkdir("SA-data");CHKERRQ(ierr);
452     }
453   }
454   ierr = PetscSNPrintf(filename,sizeof(filename),"SA-data/SA-CPS%06d.bin",id);CHKERRQ(ierr);
455   ierr = OutputBIN(filename,&viewer);CHKERRQ(ierr);
456 
457   ierr = TSGetStages(ts,&stack->numY,&Y);CHKERRQ(ierr);
458   ierr = PetscLogEventBegin(TSTrajectory_DiskWrite,ts,0,0,0);CHKERRQ(ierr);
459   ierr = WriteToDisk(stepnum,ts->ptime,ts->ptime_prev,ts->vec_sol,Y,stack->numY,stack->solution_only,viewer);CHKERRQ(ierr);
460   ierr = PetscLogEventEnd(TSTrajectory_DiskWrite,ts,0,0,0);CHKERRQ(ierr);
461   ts->trajectory->diskwrites++;
462 
463   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
464   PetscFunctionReturn(0);
465 }
466 
467 static PetscErrorCode LoadSingle(TSTrajectory tj,TS ts,Stack *stack,PetscInt id)
468 {
469   Vec            *Y;
470   PetscViewer    viewer;
471   char           filename[PETSC_MAX_PATH_LEN];
472   PetscErrorCode ierr;
473 
474   PetscFunctionBegin;
475   if (tj->monitor) {
476     ierr = PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr);
477     ierr = PetscViewerASCIIPrintf(tj->monitor,"Load a single point from file\n");CHKERRQ(ierr);
478     ierr = PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr);
479   }
480   ierr = PetscSNPrintf(filename,sizeof filename,"SA-data/SA-CPS%06d.bin",id);CHKERRQ(ierr);
481   ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);CHKERRQ(ierr);
482 
483   ierr = TSGetStages(ts,&stack->numY,&Y);CHKERRQ(ierr);
484   ierr = PetscLogEventBegin(TSTrajectory_DiskRead,ts,0,0,0);CHKERRQ(ierr);
485   ierr = ReadFromDisk(&ts->total_steps,&ts->ptime,&ts->ptime_prev,ts->vec_sol,Y,stack->numY,stack->solution_only,viewer);CHKERRQ(ierr);
486   ierr = PetscLogEventEnd(TSTrajectory_DiskRead,ts,0,0,0);CHKERRQ(ierr);
487   ts->trajectory->diskreads++;
488 
489   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
490   PetscFunctionReturn(0);
491 }
492 
493 static PetscErrorCode UpdateTS(TS ts,Stack *stack,StackElement e)
494 {
495   Vec            *Y;
496   PetscInt       i;
497   PetscErrorCode ierr;
498 
499   PetscFunctionBegin;
500   ierr = VecCopy(e->X,ts->vec_sol);CHKERRQ(ierr);
501   if (!stack->solution_only) {
502     ierr = TSGetStages(ts,&stack->numY,&Y);CHKERRQ(ierr);
503     for (i=0;i<stack->numY;i++) {
504       ierr = VecCopy(e->Y[i],Y[i]);CHKERRQ(ierr);
505     }
506   }
507   ierr = TSSetTimeStep(ts,e->timeprev-e->time);CHKERRQ(ierr); /* stepsize will be negative */
508   ts->ptime      = e->time;
509   ts->ptime_prev = e->timeprev;
510   PetscFunctionReturn(0);
511 }
512 
513 static PetscErrorCode ReCompute(TS ts,TJScheduler *tjsch,PetscInt stepnumbegin,PetscInt stepnumend)
514 {
515   Stack          *stack = &tjsch->stack;
516   PetscInt       i,adjsteps;
517   PetscErrorCode ierr;
518 
519   PetscFunctionBegin;
520   adjsteps = ts->steps;
521   ts->steps = stepnumbegin; /* global step number */
522   for (i=stepnumbegin;i<stepnumend;i++) { /* assume fixed step size */
523     if (stack->solution_only && !tjsch->skip_trajectory) { /* revolve online need this */
524       ierr = TSTrajectorySet(ts->trajectory,ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
525     }
526     ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
527     ierr = TSStep(ts);CHKERRQ(ierr);
528     if (!stack->solution_only && !tjsch->skip_trajectory) {
529       ierr = TSTrajectorySet(ts->trajectory,ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
530     }
531     ierr = TSEventHandler(ts);CHKERRQ(ierr);
532     if (!ts->steprollback) {
533       ierr = TSPostStep(ts);CHKERRQ(ierr);
534     }
535   }
536   ierr = TurnBackward(ts);CHKERRQ(ierr);
537   ts->trajectory->recomps += stepnumend-stepnumbegin; /* recomputation counter */
538   ts->steps = adjsteps;
539   ts->total_steps = stepnumend;
540   PetscFunctionReturn(0);
541 }
542 
543 static PetscErrorCode TopLevelStore(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscInt localstepnum,PetscInt laststridesize,PetscBool *done)
544 {
545   Stack          *stack = &tjsch->stack;
546   DiskStack      *diskstack = &tjsch->diskstack;
547   PetscInt       stridenum;
548   PetscErrorCode ierr;
549 
550   PetscFunctionBegin;
551   *done = PETSC_FALSE;
552   stridenum    = stepnum/tjsch->stride;
553   /* make sure saved checkpoint id starts from 1
554      skip last stride when using stridenum+1
555      skip first stride when using stridenum */
556   if (stack->solution_only) {
557     if (tjsch->save_stack) {
558       if (localstepnum == tjsch->stride-1 && stepnum < tjsch->total_steps-laststridesize) { /* current step will be saved without going through stack */
559         ierr = StackDumpAll(tj,ts,stack,stridenum+1);CHKERRQ(ierr);
560         if (tjsch->stype == TWO_LEVEL_TWO_REVOLVE) diskstack->container[++diskstack->top] = stridenum+1;
561         *done = PETSC_TRUE;
562       }
563     } else {
564       if (localstepnum == 0 && stepnum < tjsch->total_steps-laststridesize) {
565         ierr = DumpSingle(tj,ts,stack,stridenum+1);CHKERRQ(ierr);
566         if (tjsch->stype == TWO_LEVEL_TWO_REVOLVE) diskstack->container[++diskstack->top] = stridenum+1;
567         *done = PETSC_TRUE;
568       }
569     }
570   } else {
571     if (tjsch->save_stack) {
572       if (localstepnum == 0 && stepnum < tjsch->total_steps && stepnum != 0) { /* skip the first stride */
573         ierr = StackDumpAll(tj,ts,stack,stridenum);CHKERRQ(ierr);
574         if (tjsch->stype == TWO_LEVEL_TWO_REVOLVE) diskstack->container[++diskstack->top] = stridenum;
575         *done = PETSC_TRUE;
576       }
577     } else {
578       if (localstepnum == 1 && stepnum < tjsch->total_steps-laststridesize) {
579         ierr = DumpSingle(tj,ts,stack,stridenum+1);CHKERRQ(ierr);
580         if (tjsch->stype == TWO_LEVEL_TWO_REVOLVE) diskstack->container[++diskstack->top] = stridenum+1;
581         *done = PETSC_TRUE;
582       }
583     }
584   }
585   PetscFunctionReturn(0);
586 }
587 
588 static PetscErrorCode SetTrajN(TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X)
589 {
590   Stack          *stack = &tjsch->stack;
591   StackElement   e;
592   PetscErrorCode ierr;
593 
594   PetscFunctionBegin;
595   /* skip the last step */
596   if (ts->reason) { /* only affect the forward run */
597     /* update total_steps in the end of forward run */
598     if (stepnum != tjsch->total_steps) tjsch->total_steps = stepnum;
599     if (stack->solution_only) {
600       /* get rid of the solution at second last step */
601       ierr = StackPop(stack,&e);CHKERRQ(ierr);
602       ierr = ElementDestroy(stack,e);CHKERRQ(ierr);
603     }
604     PetscFunctionReturn(0);
605   }
606   /*  do not save trajectory at the recompute stage for solution_only mode */
607   if (stack->solution_only && tjsch->recompute) PetscFunctionReturn(0);
608   /* skip the first step for no_solution_only mode */
609   if (!stack->solution_only && stepnum == 0) PetscFunctionReturn(0);
610 
611   /* resize the stack */
612   if (stack->top+1 == stack->stacksize) {
613     ierr = StackResize(stack,2*stack->stacksize);CHKERRQ(ierr);
614   }
615   /* update timenext for the previous step; necessary for step adaptivity */
616   if (stack->top > -1) {
617     ierr = StackTop(stack,&e);CHKERRQ(ierr);
618     e->timenext = ts->ptime;
619   }
620   if (stepnum < stack->top) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_MEMC,"Illegal modification of a non-top stack element");
621   ierr = ElementCreate(ts,stack,&e);CHKERRQ(ierr);
622   ierr = ElementSet(ts,stack,&e,stepnum,time,X);CHKERRQ(ierr);
623   ierr = StackPush(stack,e);CHKERRQ(ierr);
624   PetscFunctionReturn(0);
625 }
626 
627 static PetscErrorCode GetTrajN(TS ts,TJScheduler *tjsch,PetscInt stepnum)
628 {
629   Stack          *stack = &tjsch->stack;
630   StackElement   e;
631   PetscErrorCode ierr;
632 
633   PetscFunctionBegin;
634   if (stepnum == tjsch->total_steps) {
635     ierr = TurnBackward(ts);CHKERRQ(ierr);
636     PetscFunctionReturn(0);
637   }
638   /* restore a checkpoint */
639   ierr = StackTop(stack,&e);CHKERRQ(ierr);
640   ierr = UpdateTS(ts,stack,e);CHKERRQ(ierr);
641   if (stack->solution_only) {/* recompute one step */
642     tjsch->recompute = PETSC_TRUE;
643     ierr = TurnForwardWithStepsize(ts,e->timenext-e->time);CHKERRQ(ierr);
644     ierr = ReCompute(ts,tjsch,e->stepnum,stepnum);CHKERRQ(ierr);
645   }
646   ierr = StackPop(stack,&e);CHKERRQ(ierr);
647   ierr = ElementDestroy(stack,e);CHKERRQ(ierr);
648   PetscFunctionReturn(0);
649 }
650 
651 static PetscErrorCode SetTrajTLNR(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X)
652 {
653   Stack          *stack = &tjsch->stack;
654   PetscInt       localstepnum,laststridesize;
655   StackElement   e;
656   PetscBool      done;
657   PetscErrorCode ierr;
658 
659   PetscFunctionBegin;
660   if (!stack->solution_only && stepnum == 0) PetscFunctionReturn(0);
661   if (stack->solution_only && stepnum == tjsch->total_steps) PetscFunctionReturn(0);
662   if (tjsch->save_stack && tjsch->recompute) PetscFunctionReturn(0);
663 
664   localstepnum = stepnum%tjsch->stride;
665   /* (stride size-1) checkpoints are saved in each stride; an extra point is added by StackDumpAll() */
666   laststridesize = tjsch->total_steps%tjsch->stride;
667   if (!laststridesize) laststridesize = tjsch->stride;
668 
669   if (!tjsch->recompute) {
670     ierr = TopLevelStore(tj,ts,tjsch,stepnum,localstepnum,laststridesize,&done);CHKERRQ(ierr);
671     if (!tjsch->save_stack && stepnum < tjsch->total_steps-laststridesize) PetscFunctionReturn(0);
672   }
673   if (!stack->solution_only && localstepnum == 0) PetscFunctionReturn(0); /* skip last point in each stride at recompute stage or last stride */
674   if (stack->solution_only && localstepnum == tjsch->stride-1) PetscFunctionReturn(0); /* skip last step in each stride at recompute stage or last stride */
675 
676   ierr = ElementCreate(ts,stack,&e);CHKERRQ(ierr);
677   ierr = ElementSet(ts,stack,&e,stepnum,time,X);CHKERRQ(ierr);
678   ierr = StackPush(stack,e);CHKERRQ(ierr);
679   PetscFunctionReturn(0);
680 }
681 
682 static PetscErrorCode GetTrajTLNR(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum)
683 {
684   Stack          *stack = &tjsch->stack;
685   PetscInt       id,localstepnum,laststridesize;
686   StackElement   e;
687   PetscErrorCode ierr;
688 
689   PetscFunctionBegin;
690   if (stepnum == tjsch->total_steps) {
691     ierr = TurnBackward(ts);CHKERRQ(ierr);
692     PetscFunctionReturn(0);
693   }
694 
695   localstepnum = stepnum%tjsch->stride;
696   laststridesize = tjsch->total_steps%tjsch->stride;
697   if (!laststridesize) laststridesize = tjsch->stride;
698   if (stack->solution_only) {
699     /* fill stack with info */
700     if (localstepnum == 0 && tjsch->total_steps-stepnum >= laststridesize) {
701       id = stepnum/tjsch->stride;
702       if (tjsch->save_stack) {
703         ierr = StackLoadAll(tj,ts,stack,id);CHKERRQ(ierr);
704         tjsch->recompute = PETSC_TRUE;
705         tjsch->skip_trajectory = PETSC_TRUE;
706         ierr = TurnForward(ts);CHKERRQ(ierr);
707         ierr = ReCompute(ts,tjsch,id*tjsch->stride-1,id*tjsch->stride);CHKERRQ(ierr);
708         tjsch->skip_trajectory = PETSC_FALSE;
709       } else {
710         ierr = LoadSingle(tj,ts,stack,id);CHKERRQ(ierr);
711         tjsch->recompute = PETSC_TRUE;
712         ierr = TurnForward(ts);CHKERRQ(ierr);
713         ierr = ReCompute(ts,tjsch,(id-1)*tjsch->stride,id*tjsch->stride);CHKERRQ(ierr);
714       }
715       PetscFunctionReturn(0);
716     }
717     /* restore a checkpoint */
718     ierr = StackPop(stack,&e);CHKERRQ(ierr);
719     ierr = UpdateTS(ts,stack,e);CHKERRQ(ierr);
720     tjsch->recompute = PETSC_TRUE;
721     tjsch->skip_trajectory = PETSC_TRUE;
722     ierr = TurnForward(ts);CHKERRQ(ierr);
723     ierr = ReCompute(ts,tjsch,e->stepnum,stepnum);CHKERRQ(ierr);
724     tjsch->skip_trajectory = PETSC_FALSE;
725     ierr = ElementDestroy(stack,e);CHKERRQ(ierr);
726   } else {
727     /* fill stack with info */
728     if (localstepnum == 0 && tjsch->total_steps-stepnum >= laststridesize) {
729       id = stepnum/tjsch->stride;
730       if (tjsch->save_stack) {
731         ierr = StackLoadAll(tj,ts,stack,id);CHKERRQ(ierr);
732       } else {
733         ierr = LoadSingle(tj,ts,stack,id);CHKERRQ(ierr);
734         ierr = ElementCreate(ts,stack,&e);CHKERRQ(ierr);
735         ierr = ElementSet(ts,stack,&e,(id-1)*tjsch->stride+1,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
736         ierr = StackPush(stack,e);CHKERRQ(ierr);
737         tjsch->recompute = PETSC_TRUE;
738         ierr = TurnForward(ts);CHKERRQ(ierr);
739         ierr = ReCompute(ts,tjsch,e->stepnum,id*tjsch->stride);CHKERRQ(ierr);
740       }
741       PetscFunctionReturn(0);
742     }
743     /* restore a checkpoint */
744     ierr = StackPop(stack,&e);CHKERRQ(ierr);
745     ierr = UpdateTS(ts,stack,e);CHKERRQ(ierr);
746     ierr = ElementDestroy(stack,e);CHKERRQ(ierr);
747   }
748   PetscFunctionReturn(0);
749 }
750 
751 #if defined(PETSC_HAVE_REVOLVE)
752 static PetscErrorCode printwhattodo(PetscViewer viewer,PetscInt whattodo,RevolveCTX *rctx,PetscInt shift)
753 {
754   PetscErrorCode ierr;
755 
756   PetscFunctionBegin;
757   if (!viewer) PetscFunctionReturn(0);
758 
759   switch(whattodo) {
760     case 1:
761       ierr = PetscViewerASCIIPrintf(viewer,"Advance from %D to %D\n",rctx->oldcapo+shift,rctx->capo+shift);CHKERRQ(ierr);
762       break;
763     case 2:
764       ierr = PetscViewerASCIIPrintf(viewer,"Store in checkpoint number %D (located in RAM)\n",rctx->check);CHKERRQ(ierr);
765       break;
766     case 3:
767       ierr = PetscViewerASCIIPrintf(viewer,"First turn: Initialize adjoints and reverse first step\n");CHKERRQ(ierr);
768       break;
769     case 4:
770       ierr = PetscViewerASCIIPrintf(viewer,"Forward and reverse one step\n");CHKERRQ(ierr);
771       break;
772     case 5:
773       ierr = PetscViewerASCIIPrintf(viewer,"Restore in checkpoint number %D (located in RAM)\n",rctx->check);CHKERRQ(ierr);
774       break;
775     case 7:
776       ierr = PetscViewerASCIIPrintf(viewer,"Store in checkpoint number %D (located on disk)\n",rctx->check);CHKERRQ(ierr);
777       break;
778     case 8:
779       ierr = PetscViewerASCIIPrintf(viewer,"Restore in checkpoint number %D (located on disk)\n",rctx->check);CHKERRQ(ierr);
780       break;
781     case -1:
782       ierr = PetscViewerASCIIPrintf(viewer,"Error!");CHKERRQ(ierr);
783       break;
784   }
785   PetscFunctionReturn(0);
786 }
787 
788 static PetscErrorCode printwhattodo2(PetscViewer viewer,PetscInt whattodo,RevolveCTX *rctx,PetscInt shift)
789 {
790   PetscErrorCode ierr;
791 
792   PetscFunctionBegin;
793   if (!viewer) PetscFunctionReturn(0);
794 
795   switch(whattodo) {
796     case 1:
797       ierr = PetscViewerASCIIPrintf(viewer,"[Top Level] Advance from stride %D to stride %D\n",rctx->oldcapo+shift,rctx->capo+shift);CHKERRQ(ierr);
798       break;
799     case 2:
800       ierr = PetscViewerASCIIPrintf(viewer,"[Top Level] Store in checkpoint number %D\n",rctx->check);CHKERRQ(ierr);
801       break;
802     case 3:
803       ierr = PetscViewerASCIIPrintf(viewer,"[Top Level] First turn: Initialize adjoints and reverse first stride\n");CHKERRQ(ierr);
804       break;
805     case 4:
806       ierr = PetscViewerASCIIPrintf(viewer,"[Top Level] Forward and reverse one stride\n");CHKERRQ(ierr);
807       break;
808     case 5:
809       ierr = PetscViewerASCIIPrintf(viewer,"[Top Level] Restore in checkpoint number %D\n",rctx->check);CHKERRQ(ierr);
810       break;
811     case 7:
812       ierr = PetscViewerASCIIPrintf(viewer,"[Top Level] Store in top-level checkpoint number %D\n",rctx->check);CHKERRQ(ierr);
813       break;
814     case 8:
815       ierr = PetscViewerASCIIPrintf(viewer,"[Top Level] Restore in top-level checkpoint number %D\n",rctx->check);CHKERRQ(ierr);
816       break;
817     case -1:
818       ierr = PetscViewerASCIIPrintf(viewer,"[Top Level] Error!");CHKERRQ(ierr);
819       break;
820   }
821   PetscFunctionReturn(0);
822 }
823 
824 static PetscErrorCode InitRevolve(PetscInt fine,PetscInt snaps,RevolveCTX *rctx)
825 {
826   PetscFunctionBegin;
827   revolve_reset();
828   revolve_create_offline(fine,snaps);
829   rctx->snaps_in       = snaps;
830   rctx->fine           = fine;
831   rctx->check          = 0;
832   rctx->capo           = 0;
833   rctx->reverseonestep = PETSC_FALSE;
834   /* check stepsleft? */
835   PetscFunctionReturn(0);
836 }
837 
838 static PetscErrorCode FastForwardRevolve(RevolveCTX *rctx)
839 {
840   PetscInt whattodo;
841 
842   PetscFunctionBegin;
843   whattodo = 0;
844   while(whattodo!=3) { /* we have to fast forward revolve to the beginning of the backward sweep due to unfriendly revolve interface */
845     whattodo = revolve_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where);
846   }
847   PetscFunctionReturn(0);
848 }
849 
850 static PetscErrorCode ApplyRevolve(PetscViewer viewer,SchedulerType stype,RevolveCTX *rctx,PetscInt total_steps,PetscInt stepnum,PetscInt localstepnum,PetscBool toplevel,PetscInt *store)
851 {
852   PetscErrorCode ierr;
853   PetscInt       shift,whattodo;
854 
855   PetscFunctionBegin;
856   *store = 0;
857   if (rctx->stepsleft > 0) { /* advance the solution without checkpointing anything as Revolve requires */
858     rctx->stepsleft--;
859     PetscFunctionReturn(0);
860   }
861   /* let Revolve determine what to do next */
862   shift         = stepnum-localstepnum;
863   rctx->oldcapo = rctx->capo;
864   rctx->capo    = localstepnum;
865 
866   if (!toplevel) whattodo = revolve_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where);
867   else whattodo = revolve2_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where);
868   if (stype == REVOLVE_ONLINE && whattodo == 8) whattodo = 5;
869   if (stype == REVOLVE_ONLINE && whattodo == 7) whattodo = 2;
870   if (!toplevel) {ierr = printwhattodo(viewer,whattodo,rctx,shift);CHKERRQ(ierr);}
871   else {ierr = printwhattodo2(viewer,whattodo,rctx,shift);CHKERRQ(ierr);}
872   if (whattodo == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in the Revolve library");
873   if (whattodo == 1) { /* advance some time steps */
874     if (stype == REVOLVE_ONLINE && rctx->capo >= total_steps-1) {
875       revolve_turn(total_steps,&rctx->capo,&rctx->fine);
876       if (!toplevel) {ierr = printwhattodo(viewer,whattodo,rctx,shift);CHKERRQ(ierr);}
877       else {ierr = printwhattodo2(viewer,whattodo,rctx,shift);CHKERRQ(ierr);}
878     }
879     rctx->stepsleft = rctx->capo-rctx->oldcapo-1;
880   }
881   if (whattodo == 3 || whattodo == 4) { /* ready for a reverse step */
882     rctx->reverseonestep = PETSC_TRUE;
883   }
884   if (whattodo == 5) { /* restore a checkpoint and ask Revolve what to do next */
885     rctx->oldcapo = rctx->capo;
886     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*/
887     else whattodo = revolve2_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where);
888     if (!toplevel) {ierr = printwhattodo(viewer,whattodo,rctx,shift);CHKERRQ(ierr);}
889     else {ierr = printwhattodo2(viewer,whattodo,rctx,shift);CHKERRQ(ierr);}
890     if (whattodo == 3 || whattodo == 4) rctx->reverseonestep = PETSC_TRUE;
891     if (whattodo == 1) rctx->stepsleft = rctx->capo-rctx->oldcapo;
892   }
893   if (whattodo == 7) { /* save the checkpoint to disk */
894     *store = 2;
895     rctx->oldcapo = rctx->capo;
896     whattodo = revolve_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where); /* must return 1 */
897     ierr = printwhattodo(viewer,whattodo,rctx,shift);CHKERRQ(ierr);
898     rctx->stepsleft = rctx->capo-rctx->oldcapo-1;
899   }
900   if (whattodo == 2) { /* store a checkpoint to RAM and ask Revolve how many time steps to advance next */
901     *store = 1;
902     rctx->oldcapo = rctx->capo;
903     if (!toplevel) whattodo = revolve_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where); /* must return 1 */
904     else whattodo = revolve2_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where);
905     if (!toplevel) {ierr = printwhattodo(viewer,whattodo,rctx,shift);CHKERRQ(ierr);}
906     else {ierr = printwhattodo2(viewer,whattodo,rctx,shift);CHKERRQ(ierr);}
907     if (stype == REVOLVE_ONLINE && rctx->capo >= total_steps-1) {
908       revolve_turn(total_steps,&rctx->capo,&rctx->fine);
909       ierr = printwhattodo(viewer,whattodo,rctx,shift);CHKERRQ(ierr);
910     }
911     rctx->stepsleft = rctx->capo-rctx->oldcapo-1;
912   }
913   PetscFunctionReturn(0);
914 }
915 
916 static PetscErrorCode SetTrajROF(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X)
917 {
918   Stack          *stack = &tjsch->stack;
919   PetscInt       store;
920   StackElement   e;
921   PetscErrorCode ierr;
922 
923   PetscFunctionBegin;
924   if (!stack->solution_only && stepnum == 0) PetscFunctionReturn(0);
925   if (stack->solution_only && stepnum == tjsch->total_steps) PetscFunctionReturn(0);
926   ierr = ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,tjsch->total_steps,stepnum,stepnum,PETSC_FALSE,&store);CHKERRQ(ierr);
927   if (store == 1) {
928     if (stepnum < stack->top) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_MEMC,"Illegal modification of a non-top stack element");
929     ierr = ElementCreate(ts,stack,&e);CHKERRQ(ierr);
930     ierr = ElementSet(ts,stack,&e,stepnum,time,X);CHKERRQ(ierr);
931     ierr = StackPush(stack,e);CHKERRQ(ierr);
932   }
933   PetscFunctionReturn(0);
934 }
935 
936 static PetscErrorCode GetTrajROF(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum)
937 {
938   Stack          *stack = &tjsch->stack;
939   PetscInt       whattodo,shift,store;
940   StackElement   e;
941   PetscErrorCode ierr;
942 
943   PetscFunctionBegin;
944   if (stepnum == 0 || stepnum == tjsch->total_steps) {
945     ierr = TurnBackward(ts);CHKERRQ(ierr);
946     tjsch->rctx->reverseonestep = PETSC_FALSE;
947     PetscFunctionReturn(0);
948   }
949   /* restore a checkpoint */
950   ierr = StackTop(stack,&e);CHKERRQ(ierr);
951   ierr = UpdateTS(ts,stack,e);CHKERRQ(ierr);
952   if (stack->solution_only) { /* start with restoring a checkpoint */
953     tjsch->rctx->capo = stepnum;
954     tjsch->rctx->oldcapo = tjsch->rctx->capo;
955     shift = 0;
956     whattodo = revolve_action(&tjsch->rctx->check,&tjsch->rctx->capo,&tjsch->rctx->fine,tjsch->rctx->snaps_in,&tjsch->rctx->info,&tjsch->rctx->where);
957     ierr = printwhattodo(tj->monitor,whattodo,tjsch->rctx,shift);CHKERRQ(ierr);
958   } else { /* 2 revolve actions: restore a checkpoint and then advance */
959     ierr = ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,tjsch->total_steps,stepnum,stepnum,PETSC_FALSE,&store);CHKERRQ(ierr);
960     if (tj->monitor) {
961       ierr = PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr);
962       ierr = PetscViewerASCIIPrintf(tj->monitor,"Skip the step from %D to %D (stage values already checkpointed)\n",tjsch->rctx->oldcapo,tjsch->rctx->oldcapo+1);CHKERRQ(ierr);
963       ierr = PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr);
964     }
965     if (!tjsch->rctx->reverseonestep && tjsch->rctx->stepsleft > 0) tjsch->rctx->stepsleft--;
966   }
967   if (stack->solution_only || (!stack->solution_only && e->stepnum < stepnum)) {
968     tjsch->recompute = PETSC_TRUE;
969     ierr = TurnForward(ts);CHKERRQ(ierr);
970     ierr = ReCompute(ts,tjsch,e->stepnum,stepnum);CHKERRQ(ierr);
971   }
972   if ((stack->solution_only && e->stepnum+1 == stepnum) || (!stack->solution_only && e->stepnum == stepnum)) {
973     ierr = StackPop(stack,&e);CHKERRQ(ierr);
974     ierr = ElementDestroy(stack,e);CHKERRQ(ierr);
975   }
976   tjsch->rctx->reverseonestep = PETSC_FALSE;
977   PetscFunctionReturn(0);
978 }
979 
980 static PetscErrorCode SetTrajRON(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X)
981 {
982   Stack          *stack = &tjsch->stack;
983   Vec            *Y;
984   PetscInt       i,store;
985   PetscReal      timeprev;
986   StackElement   e;
987   RevolveCTX     *rctx = tjsch->rctx;
988   PetscErrorCode ierr;
989 
990   PetscFunctionBegin;
991   if (!stack->solution_only && stepnum == 0) PetscFunctionReturn(0);
992   if (stack->solution_only && stepnum == tjsch->total_steps) PetscFunctionReturn(0);
993   ierr = ApplyRevolve(tj->monitor,tjsch->stype,rctx,tjsch->total_steps,stepnum,stepnum,PETSC_FALSE,&store);CHKERRQ(ierr);
994   if (store == 1) {
995     if (rctx->check != stack->top+1) { /* overwrite some non-top checkpoint in the stack */
996       ierr = StackFind(stack,&e,rctx->check);CHKERRQ(ierr);
997       ierr = VecCopy(X,e->X);CHKERRQ(ierr);
998       if (stack->numY > 0 && !stack->solution_only) {
999         ierr = TSGetStages(ts,&stack->numY,&Y);CHKERRQ(ierr);
1000         for (i=0;i<stack->numY;i++) {
1001           ierr = VecCopy(Y[i],e->Y[i]);CHKERRQ(ierr);
1002         }
1003       }
1004       e->stepnum  = stepnum;
1005       e->time     = time;
1006       ierr        = TSGetPrevTime(ts,&timeprev);CHKERRQ(ierr);
1007       e->timeprev = timeprev;
1008     } else {
1009       if (stepnum < stack->top) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_MEMC,"Illegal modification of a non-top stack element");
1010       ierr = ElementCreate(ts,stack,&e);CHKERRQ(ierr);
1011       ierr = ElementSet(ts,stack,&e,stepnum,time,X);CHKERRQ(ierr);
1012       ierr = StackPush(stack,e);CHKERRQ(ierr);
1013     }
1014   }
1015   PetscFunctionReturn(0);
1016 }
1017 
1018 static PetscErrorCode GetTrajRON(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum)
1019 {
1020   Stack          *stack = &tjsch->stack;
1021   PetscInt       whattodo,shift;
1022   StackElement   e;
1023   PetscErrorCode ierr;
1024 
1025   PetscFunctionBegin;
1026   if (stepnum == 0 || stepnum == tjsch->total_steps) {
1027     ierr = TurnBackward(ts);CHKERRQ(ierr);
1028     tjsch->rctx->reverseonestep = PETSC_FALSE;
1029     PetscFunctionReturn(0);
1030   }
1031   tjsch->rctx->capo = stepnum;
1032   tjsch->rctx->oldcapo = tjsch->rctx->capo;
1033   shift = 0;
1034   whattodo = revolve_action(&tjsch->rctx->check,&tjsch->rctx->capo,&tjsch->rctx->fine,tjsch->rctx->snaps_in,&tjsch->rctx->info,&tjsch->rctx->where); /* whattodo=restore */
1035   if (whattodo == 8) whattodo = 5;
1036   ierr = printwhattodo(tj->monitor,whattodo,tjsch->rctx,shift);CHKERRQ(ierr);
1037   /* restore a checkpoint */
1038   ierr = StackFind(stack,&e,tjsch->rctx->check);CHKERRQ(ierr);
1039   ierr = UpdateTS(ts,stack,e);CHKERRQ(ierr);
1040   if (!stack->solution_only) { /* whattodo must be 5 */
1041     /* ask Revolve what to do next */
1042     tjsch->rctx->oldcapo = tjsch->rctx->capo;
1043     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*/
1044     ierr = printwhattodo(tj->monitor,whattodo,tjsch->rctx,shift);CHKERRQ(ierr);
1045     if (whattodo == 3 || whattodo == 4) tjsch->rctx->reverseonestep = PETSC_TRUE;
1046     if (whattodo == 1) tjsch->rctx->stepsleft = tjsch->rctx->capo-tjsch->rctx->oldcapo;
1047     if (tj->monitor) {
1048       ierr = PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr);
1049       ierr = PetscViewerASCIIPrintf(tj->monitor,"Skip the step from %D to %D (stage values already checkpointed)\n",tjsch->rctx->oldcapo,tjsch->rctx->oldcapo+1);CHKERRQ(ierr);
1050       ierr = PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr);
1051     }
1052     if (!tjsch->rctx->reverseonestep && tjsch->rctx->stepsleft > 0) tjsch->rctx->stepsleft--;
1053   }
1054   if (stack->solution_only || (!stack->solution_only && e->stepnum < stepnum)) {
1055     tjsch->recompute = PETSC_TRUE;
1056     ierr = TurnForward(ts);CHKERRQ(ierr);
1057     ierr = ReCompute(ts,tjsch,e->stepnum,stepnum);CHKERRQ(ierr);
1058   }
1059   tjsch->rctx->reverseonestep = PETSC_FALSE;
1060   PetscFunctionReturn(0);
1061 }
1062 
1063 static PetscErrorCode SetTrajTLR(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X)
1064 {
1065   Stack          *stack = &tjsch->stack;
1066   PetscInt       store,localstepnum,laststridesize;
1067   StackElement   e;
1068   PetscBool      done = PETSC_FALSE;
1069   PetscErrorCode ierr;
1070 
1071   PetscFunctionBegin;
1072   if (!stack->solution_only && stepnum == 0) PetscFunctionReturn(0);
1073   if (stack->solution_only && stepnum == tjsch->total_steps) PetscFunctionReturn(0);
1074 
1075   localstepnum = stepnum%tjsch->stride;
1076   laststridesize = tjsch->total_steps%tjsch->stride;
1077   if (!laststridesize) laststridesize = tjsch->stride;
1078 
1079   if (!tjsch->recompute) {
1080     ierr = TopLevelStore(tj,ts,tjsch,stepnum,localstepnum,laststridesize,&done);CHKERRQ(ierr);
1081     /* revolve is needed for the last stride; different starting points for last stride between solutin_only and !solutin_only */
1082     if (!stack->solution_only && !tjsch->save_stack && stepnum <= tjsch->total_steps-laststridesize) PetscFunctionReturn(0);
1083     if (stack->solution_only && !tjsch->save_stack && stepnum < tjsch->total_steps-laststridesize) PetscFunctionReturn(0);
1084   }
1085   if (tjsch->save_stack && done) {
1086     ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr);
1087     PetscFunctionReturn(0);
1088   }
1089   if (laststridesize < tjsch->stride) {
1090     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 */
1091       ierr = InitRevolve(laststridesize,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr);
1092     }
1093     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 */
1094       ierr = InitRevolve(laststridesize,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr);
1095     }
1096   }
1097   ierr = ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,tjsch->total_steps,stepnum,localstepnum,PETSC_FALSE,&store);CHKERRQ(ierr);
1098   if (store == 1) {
1099     if (localstepnum < stack->top) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_MEMC,"Illegal modification of a non-top stack element");
1100     ierr = ElementCreate(ts,stack,&e);CHKERRQ(ierr);
1101     ierr = ElementSet(ts,stack,&e,stepnum,time,X);CHKERRQ(ierr);
1102     ierr = StackPush(stack,e);CHKERRQ(ierr);
1103   }
1104   PetscFunctionReturn(0);
1105 }
1106 
1107 static PetscErrorCode GetTrajTLR(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum)
1108 {
1109   Stack          *stack = &tjsch->stack;
1110   PetscInt       whattodo,shift;
1111   PetscInt       localstepnum,stridenum,laststridesize,store;
1112   StackElement   e;
1113   PetscErrorCode ierr;
1114 
1115   PetscFunctionBegin;
1116   localstepnum = stepnum%tjsch->stride;
1117   stridenum    = stepnum/tjsch->stride;
1118   if (stepnum == tjsch->total_steps) {
1119     ierr = TurnBackward(ts);CHKERRQ(ierr);
1120     tjsch->rctx->reverseonestep = PETSC_FALSE;
1121     PetscFunctionReturn(0);
1122   }
1123   laststridesize = tjsch->total_steps%tjsch->stride;
1124   if (!laststridesize) laststridesize = tjsch->stride;
1125   if (stack->solution_only) {
1126     /* fill stack */
1127     if (localstepnum == 0 && stepnum <= tjsch->total_steps-laststridesize) {
1128       if (tjsch->save_stack) {
1129         ierr = StackLoadAll(tj,ts,stack,stridenum);CHKERRQ(ierr);
1130         ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr);
1131         ierr = FastForwardRevolve(tjsch->rctx);CHKERRQ(ierr);
1132         tjsch->recompute = PETSC_TRUE;
1133         tjsch->skip_trajectory = PETSC_TRUE;
1134         ierr = TurnForward(ts);CHKERRQ(ierr);
1135         ierr = ReCompute(ts,tjsch,stridenum*tjsch->stride-1,stridenum*tjsch->stride);CHKERRQ(ierr);
1136         tjsch->skip_trajectory = PETSC_FALSE;
1137       } else {
1138         ierr = LoadSingle(tj,ts,stack,stridenum);CHKERRQ(ierr);
1139         ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr);
1140         tjsch->recompute = PETSC_TRUE;
1141         ierr = TurnForward(ts);CHKERRQ(ierr);
1142         ierr = ReCompute(ts,tjsch,(stridenum-1)*tjsch->stride,stridenum*tjsch->stride);CHKERRQ(ierr);
1143       }
1144       PetscFunctionReturn(0);
1145     }
1146     /* restore a checkpoint */
1147     ierr = StackTop(stack,&e);CHKERRQ(ierr);
1148     ierr = UpdateTS(ts,stack,e);CHKERRQ(ierr);
1149     /* start with restoring a checkpoint */
1150     tjsch->rctx->capo = stepnum;
1151     tjsch->rctx->oldcapo = tjsch->rctx->capo;
1152     shift = stepnum-localstepnum;
1153     whattodo = revolve_action(&tjsch->rctx->check,&tjsch->rctx->capo,&tjsch->rctx->fine,tjsch->rctx->snaps_in,&tjsch->rctx->info,&tjsch->rctx->where);
1154     ierr = printwhattodo(tj->monitor,whattodo,tjsch->rctx,shift);CHKERRQ(ierr);
1155     tjsch->recompute = PETSC_TRUE;
1156     ierr = TurnForward(ts);CHKERRQ(ierr);
1157     ierr = ReCompute(ts,tjsch,e->stepnum,stepnum);CHKERRQ(ierr);
1158     if (e->stepnum+1 == stepnum) {
1159       ierr = StackPop(stack,&e);CHKERRQ(ierr);
1160       ierr = ElementDestroy(stack,e);CHKERRQ(ierr);
1161     }
1162   } else {
1163     /* fill stack with info */
1164     if (localstepnum == 0 && tjsch->total_steps-stepnum >= laststridesize) {
1165       if (tjsch->save_stack) {
1166         ierr = StackLoadAll(tj,ts,stack,stridenum);CHKERRQ(ierr);
1167         ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr);
1168         ierr = FastForwardRevolve(tjsch->rctx);CHKERRQ(ierr);
1169       } else {
1170         ierr = LoadSingle(tj,ts,stack,stridenum);CHKERRQ(ierr);
1171         ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr);
1172         ierr = ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,tjsch->total_steps,(stridenum-1)*tjsch->stride+1,1,PETSC_FALSE,&store);CHKERRQ(ierr);
1173         if (tj->monitor) {
1174           ierr = PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr);
1175           ierr = PetscViewerASCIIPrintf(tj->monitor,"Skip the step from %D to %D (stage values already checkpointed)\n",(stridenum-1)*tjsch->stride+tjsch->rctx->oldcapo,(stridenum-1)*tjsch->stride+tjsch->rctx->oldcapo+1);CHKERRQ(ierr);
1176           ierr = PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr);
1177         }
1178         ierr = ElementCreate(ts,stack,&e);CHKERRQ(ierr);
1179         ierr = ElementSet(ts,stack,&e,(stridenum-1)*tjsch->stride+1,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
1180         ierr = StackPush(stack,e);CHKERRQ(ierr);
1181         tjsch->recompute = PETSC_TRUE;
1182         ierr = TurnForward(ts);CHKERRQ(ierr);
1183         ierr = ReCompute(ts,tjsch,e->stepnum,stridenum*tjsch->stride);CHKERRQ(ierr);
1184       }
1185       PetscFunctionReturn(0);
1186     }
1187     /* restore a checkpoint */
1188     ierr = StackTop(stack,&e);CHKERRQ(ierr);
1189     ierr = UpdateTS(ts,stack,e);CHKERRQ(ierr);
1190     /* 2 revolve actions: restore a checkpoint and then advance */
1191     ierr = ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,tjsch->total_steps,stepnum,localstepnum,PETSC_FALSE,&store);CHKERRQ(ierr);
1192     if (tj->monitor) {
1193       ierr = PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr);
1194       ierr = PetscViewerASCIIPrintf(tj->monitor,"Skip the step from %D to %D (stage values already checkpointed)\n",stepnum-localstepnum+tjsch->rctx->oldcapo,stepnum-localstepnum+tjsch->rctx->oldcapo+1);CHKERRQ(ierr);
1195       ierr = PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr);
1196     }
1197     if (!tjsch->rctx->reverseonestep && tjsch->rctx->stepsleft > 0) tjsch->rctx->stepsleft--;
1198     if (e->stepnum < stepnum) {
1199       tjsch->recompute = PETSC_TRUE;
1200       ierr = TurnForward(ts);CHKERRQ(ierr);
1201       ierr = ReCompute(ts,tjsch,e->stepnum,stepnum);CHKERRQ(ierr);
1202     }
1203     if (e->stepnum == stepnum) {
1204       ierr = StackPop(stack,&e);CHKERRQ(ierr);
1205       ierr = ElementDestroy(stack,e);CHKERRQ(ierr);
1206     }
1207   }
1208   tjsch->rctx->reverseonestep = PETSC_FALSE;
1209   PetscFunctionReturn(0);
1210 }
1211 
1212 static PetscErrorCode SetTrajTLTR(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X)
1213 {
1214   Stack          *stack = &tjsch->stack;
1215   PetscInt       store,localstepnum,stridenum,laststridesize;
1216   StackElement   e;
1217   PetscBool      done = PETSC_FALSE;
1218   PetscErrorCode ierr;
1219 
1220   PetscFunctionBegin;
1221   if (!stack->solution_only && stepnum == 0) PetscFunctionReturn(0);
1222   if (stack->solution_only && stepnum == tjsch->total_steps) PetscFunctionReturn(0);
1223 
1224   localstepnum = stepnum%tjsch->stride; /* index at the bottom level (inside a stride) */
1225   stridenum    = stepnum/tjsch->stride; /* index at the top level */
1226   laststridesize = tjsch->total_steps%tjsch->stride;
1227   if (!laststridesize) laststridesize = tjsch->stride;
1228   if (stack->solution_only && localstepnum == 0 && !tjsch->rctx2->reverseonestep) {
1229     ierr = ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx2,(tjsch->total_steps+tjsch->stride-1)/tjsch->stride,stridenum,stridenum,PETSC_TRUE,&tjsch->store_stride);CHKERRQ(ierr);
1230     if (laststridesize < tjsch->stride && stepnum == tjsch->total_steps-laststridesize) {
1231       ierr = InitRevolve(laststridesize,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr);
1232     }
1233   }
1234   if (!stack->solution_only && localstepnum == 1 && !tjsch->rctx2->reverseonestep) {
1235     ierr = ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx2,(tjsch->total_steps+tjsch->stride-1)/tjsch->stride,stridenum,stridenum,PETSC_TRUE,&tjsch->store_stride);CHKERRQ(ierr);
1236     if (laststridesize < tjsch->stride && stepnum == tjsch->total_steps-laststridesize+1) {
1237       ierr = InitRevolve(laststridesize,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr);
1238     }
1239   }
1240   if (tjsch->store_stride) {
1241     ierr = TopLevelStore(tj,ts,tjsch,stepnum,localstepnum,laststridesize,&done);CHKERRQ(ierr);
1242     if (done) {
1243       ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr);
1244       PetscFunctionReturn(0);
1245     }
1246   }
1247   if (stepnum < tjsch->total_steps-laststridesize) {
1248     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 */
1249     if (!tjsch->save_stack && !tjsch->rctx2->reverseonestep) PetscFunctionReturn(0); /* store operation does not require revolve be called at bottom level */
1250   }
1251   /* 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 */
1252   if (!stack->solution_only && localstepnum == 0 && stepnum != tjsch->total_steps && !tjsch->recompute) PetscFunctionReturn(0);
1253   ierr = ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,tjsch->total_steps,stepnum,localstepnum,PETSC_FALSE,&store);CHKERRQ(ierr);
1254   if (store == 1) {
1255     if (localstepnum < stack->top) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_MEMC,"Illegal modification of a non-top stack element");
1256     ierr = ElementCreate(ts,stack,&e);CHKERRQ(ierr);
1257     ierr = ElementSet(ts,stack,&e,stepnum,time,X);CHKERRQ(ierr);
1258     ierr = StackPush(stack,e);CHKERRQ(ierr);
1259   }
1260   PetscFunctionReturn(0);
1261 }
1262 
1263 static PetscErrorCode GetTrajTLTR(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum)
1264 {
1265   Stack          *stack = &tjsch->stack;
1266   DiskStack      *diskstack = &tjsch->diskstack;
1267   PetscInt       whattodo,shift;
1268   PetscInt       localstepnum,stridenum,restoredstridenum,laststridesize,store;
1269   StackElement   e;
1270   PetscErrorCode ierr;
1271 
1272   PetscFunctionBegin;
1273   localstepnum = stepnum%tjsch->stride;
1274   stridenum    = stepnum/tjsch->stride;
1275   if (stepnum == tjsch->total_steps) {
1276     ierr = TurnBackward(ts);CHKERRQ(ierr);
1277     tjsch->rctx->reverseonestep = PETSC_FALSE;
1278     PetscFunctionReturn(0);
1279   }
1280   laststridesize = tjsch->total_steps%tjsch->stride;
1281   if (!laststridesize) laststridesize = tjsch->stride;
1282   /*
1283    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:
1284      Case 1 (save_stack)
1285        Restore a disk checkpoint; update TS with the last element in the restored data; recompute to the current point.
1286      Case 2 (!save_stack)
1287        Restore a disk checkpoint; update TS with the restored point; recompute to the current point.
1288   */
1289   if (localstepnum == 0 && stepnum <= tjsch->total_steps-laststridesize) {
1290     /* restore the top element in the stack for disk checkpoints */
1291     restoredstridenum = diskstack->container[diskstack->top];
1292     tjsch->rctx2->reverseonestep = PETSC_FALSE;
1293     /* top-level revolve must be applied before current step, just like the solution_only mode for single-level revolve */
1294     if (!tjsch->save_stack && stack->solution_only) { /* start with restoring a checkpoint */
1295       tjsch->rctx2->capo = stridenum;
1296       tjsch->rctx2->oldcapo = tjsch->rctx2->capo;
1297       shift = 0;
1298       whattodo = revolve2_action(&tjsch->rctx2->check,&tjsch->rctx2->capo,&tjsch->rctx2->fine,tjsch->rctx2->snaps_in,&tjsch->rctx2->info,&tjsch->rctx2->where);
1299       ierr = printwhattodo2(tj->monitor,whattodo,tjsch->rctx2,shift);CHKERRQ(ierr);
1300     } else { /* 2 revolve actions: restore a checkpoint and then advance */
1301       ierr = ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx2,(tjsch->total_steps+tjsch->stride-1)/tjsch->stride,stridenum,stridenum,PETSC_TRUE,&tjsch->store_stride);CHKERRQ(ierr);
1302       if (tj->monitor) {
1303         ierr = PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr);
1304         ierr = PetscViewerASCIIPrintf(tj->monitor,"[Top Level] Skip the stride from %D to %D (stage values already checkpointed)\n",tjsch->rctx2->oldcapo,tjsch->rctx2->oldcapo+1);CHKERRQ(ierr);
1305         ierr = PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr);
1306       }
1307       if (!tjsch->rctx2->reverseonestep && tjsch->rctx2->stepsleft > 0) tjsch->rctx2->stepsleft--;
1308     }
1309     /* fill stack */
1310     if (stack->solution_only) {
1311       if (tjsch->save_stack) {
1312         if (restoredstridenum < stridenum) {
1313           ierr = StackLoadLast(tj,ts,stack,restoredstridenum);CHKERRQ(ierr);
1314         } else {
1315           ierr = StackLoadAll(tj,ts,stack,restoredstridenum);CHKERRQ(ierr);
1316         }
1317         /* recompute one step ahead */
1318         tjsch->recompute = PETSC_TRUE;
1319         tjsch->skip_trajectory = PETSC_TRUE;
1320         ierr = TurnForward(ts);CHKERRQ(ierr);
1321         ierr = ReCompute(ts,tjsch,stridenum*tjsch->stride-1,stridenum*tjsch->stride);CHKERRQ(ierr);
1322         tjsch->skip_trajectory = PETSC_FALSE;
1323         if (restoredstridenum < stridenum) {
1324           ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr);
1325           tjsch->recompute = PETSC_TRUE;
1326           ierr = TurnForward(ts);CHKERRQ(ierr);
1327           ierr = ReCompute(ts,tjsch,restoredstridenum*tjsch->stride,stepnum);CHKERRQ(ierr);
1328         } else { /* stack ready, fast forward revolve status */
1329           ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr);
1330           ierr = FastForwardRevolve(tjsch->rctx);CHKERRQ(ierr);
1331         }
1332       } else {
1333         ierr = LoadSingle(tj,ts,stack,restoredstridenum);CHKERRQ(ierr);
1334         ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr);
1335         tjsch->recompute = PETSC_TRUE;
1336         ierr = TurnForward(ts);CHKERRQ(ierr);
1337         ierr = ReCompute(ts,tjsch,(restoredstridenum-1)*tjsch->stride,stepnum);CHKERRQ(ierr);
1338       }
1339     } else {
1340       if (tjsch->save_stack) {
1341         if (restoredstridenum < stridenum) {
1342           ierr = StackLoadLast(tj,ts,stack,restoredstridenum);CHKERRQ(ierr);
1343           /* reset revolve */
1344           ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr);
1345           tjsch->recompute = PETSC_TRUE;
1346           ierr = TurnForward(ts);CHKERRQ(ierr);
1347           ierr = ReCompute(ts,tjsch,restoredstridenum*tjsch->stride,stepnum);CHKERRQ(ierr);
1348         } else { /* stack ready, fast forward revolve status */
1349           ierr = StackLoadAll(tj,ts,stack,restoredstridenum);CHKERRQ(ierr);
1350           ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr);
1351           ierr = FastForwardRevolve(tjsch->rctx);CHKERRQ(ierr);
1352         }
1353       } else {
1354         ierr = LoadSingle(tj,ts,stack,restoredstridenum);CHKERRQ(ierr);
1355         ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr);
1356         /* push first element to stack */
1357         if (tjsch->store_stride || tjsch->rctx2->reverseonestep) {
1358           shift = (restoredstridenum-1)*tjsch->stride-localstepnum;
1359           ierr = ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,tjsch->total_steps,(restoredstridenum-1)*tjsch->stride+1,1,PETSC_FALSE,&store);CHKERRQ(ierr);
1360           if (tj->monitor) {
1361             ierr = PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr);
1362             ierr = PetscViewerASCIIPrintf(tj->monitor,"Skip the step from %D to %D (stage values already checkpointed)\n",(restoredstridenum-1)*tjsch->stride,(restoredstridenum-1)*tjsch->stride+1);CHKERRQ(ierr);
1363             ierr = PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr);
1364           }
1365           ierr = ElementCreate(ts,stack,&e);CHKERRQ(ierr);
1366           ierr = ElementSet(ts,stack,&e,(restoredstridenum-1)*tjsch->stride+1,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
1367           ierr = StackPush(stack,e);CHKERRQ(ierr);
1368         }
1369         tjsch->recompute = PETSC_TRUE;
1370         ierr = TurnForward(ts);CHKERRQ(ierr);
1371         ierr = ReCompute(ts,tjsch,(restoredstridenum-1)*tjsch->stride+1,stepnum);CHKERRQ(ierr);
1372       }
1373     }
1374     if (restoredstridenum == stridenum) diskstack->top--;
1375     tjsch->rctx->reverseonestep = PETSC_FALSE;
1376     PetscFunctionReturn(0);
1377   }
1378 
1379   if (stack->solution_only) {
1380     /* restore a checkpoint */
1381     ierr = StackTop(stack,&e);CHKERRQ(ierr);
1382     ierr = UpdateTS(ts,stack,e);CHKERRQ(ierr);
1383     /* start with restoring a checkpoint */
1384     tjsch->rctx->capo = stepnum;
1385     tjsch->rctx->oldcapo = tjsch->rctx->capo;
1386     shift = stepnum-localstepnum;
1387     whattodo = revolve_action(&tjsch->rctx->check,&tjsch->rctx->capo,&tjsch->rctx->fine,tjsch->rctx->snaps_in,&tjsch->rctx->info,&tjsch->rctx->where);
1388     ierr = printwhattodo(tj->monitor,whattodo,tjsch->rctx,shift);CHKERRQ(ierr);
1389     tjsch->recompute = PETSC_TRUE;
1390     ierr = TurnForward(ts);CHKERRQ(ierr);
1391     ierr = ReCompute(ts,tjsch,e->stepnum,stepnum);CHKERRQ(ierr);
1392     if (e->stepnum+1 == stepnum) {
1393       ierr = StackPop(stack,&e);CHKERRQ(ierr);
1394       ierr = ElementDestroy(stack,e);CHKERRQ(ierr);
1395     }
1396   } else {
1397     /* restore a checkpoint */
1398     ierr = StackTop(stack,&e);CHKERRQ(ierr);
1399     ierr = UpdateTS(ts,stack,e);CHKERRQ(ierr);
1400     /* 2 revolve actions: restore a checkpoint and then advance */
1401     ierr = ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,tjsch->total_steps,stepnum,localstepnum,PETSC_FALSE,&store);CHKERRQ(ierr);
1402     if (tj->monitor) {
1403       ierr = PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr);
1404       ierr = PetscViewerASCIIPrintf(tj->monitor,"Skip the step from %D to %D (stage values already checkpointed)\n",stepnum-localstepnum+tjsch->rctx->oldcapo,stepnum-localstepnum+tjsch->rctx->oldcapo+1);CHKERRQ(ierr);
1405       ierr = PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr);
1406     }
1407     if (!tjsch->rctx->reverseonestep && tjsch->rctx->stepsleft > 0) tjsch->rctx->stepsleft--;
1408     if (e->stepnum < stepnum) {
1409       tjsch->recompute = PETSC_TRUE;
1410       ierr = TurnForward(ts);CHKERRQ(ierr);
1411       ierr = ReCompute(ts,tjsch,e->stepnum,stepnum);CHKERRQ(ierr);
1412     }
1413     if (e->stepnum == stepnum) {
1414       ierr = StackPop(stack,&e);CHKERRQ(ierr);
1415       ierr = ElementDestroy(stack,e);CHKERRQ(ierr);
1416     }
1417   }
1418   tjsch->rctx->reverseonestep = PETSC_FALSE;
1419   PetscFunctionReturn(0);
1420 }
1421 
1422 static PetscErrorCode SetTrajRMS(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X)
1423 {
1424   Stack          *stack = &tjsch->stack;
1425   PetscInt       store;
1426   StackElement   e;
1427   PetscErrorCode ierr;
1428 
1429   PetscFunctionBegin;
1430   if (!stack->solution_only && stepnum == 0) PetscFunctionReturn(0);
1431   if (stack->solution_only && stepnum == tjsch->total_steps) PetscFunctionReturn(0);
1432   ierr = ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,tjsch->total_steps,stepnum,stepnum,PETSC_FALSE,&store);CHKERRQ(ierr);
1433   if (store == 1){
1434     if (stepnum < stack->top) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_MEMC,"Illegal modification of a non-top stack element");
1435     ierr = ElementCreate(ts,stack,&e);CHKERRQ(ierr);
1436     ierr = ElementSet(ts,stack,&e,stepnum,time,X);CHKERRQ(ierr);
1437     ierr = StackPush(stack,e);CHKERRQ(ierr);
1438   } else if (store == 2) {
1439     ierr = DumpSingle(tj,ts,stack,tjsch->rctx->check+1);CHKERRQ(ierr);
1440   }
1441   PetscFunctionReturn(0);
1442 }
1443 
1444 static PetscErrorCode GetTrajRMS(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum)
1445 {
1446   Stack          *stack = &tjsch->stack;
1447   PetscInt       whattodo,shift;
1448   PetscInt       restart;
1449   PetscBool      ondisk;
1450   StackElement   e;
1451   PetscErrorCode ierr;
1452 
1453   PetscFunctionBegin;
1454   if (stepnum == 0 || stepnum == tjsch->total_steps) {
1455     ierr = TurnBackward(ts);CHKERRQ(ierr);
1456     tjsch->rctx->reverseonestep = PETSC_FALSE;
1457     PetscFunctionReturn(0);
1458   }
1459   tjsch->rctx->capo = stepnum;
1460   tjsch->rctx->oldcapo = tjsch->rctx->capo;
1461   shift = 0;
1462   whattodo = revolve_action(&tjsch->rctx->check,&tjsch->rctx->capo,&tjsch->rctx->fine,tjsch->rctx->snaps_in,&tjsch->rctx->info,&tjsch->rctx->where); /* whattodo=restore */
1463   ierr = printwhattodo(tj->monitor,whattodo,tjsch->rctx,shift);CHKERRQ(ierr);
1464   /* restore a checkpoint */
1465   restart = tjsch->rctx->capo;
1466   if (!tjsch->rctx->where) {
1467     ondisk = PETSC_TRUE;
1468     ierr = LoadSingle(tj,ts,stack,tjsch->rctx->check+1);CHKERRQ(ierr);
1469     ierr = TurnBackward(ts);CHKERRQ(ierr);
1470   } else {
1471     ondisk = PETSC_FALSE;
1472     ierr = StackTop(stack,&e);CHKERRQ(ierr);
1473     ierr = UpdateTS(ts,stack,e);CHKERRQ(ierr);
1474   }
1475   if (!stack->solution_only) { /* whattodo must be 5 or 8 */
1476     /* ask Revolve what to do next */
1477     tjsch->rctx->oldcapo = tjsch->rctx->capo;
1478     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*/
1479     ierr = printwhattodo(tj->monitor,whattodo,tjsch->rctx,shift);CHKERRQ(ierr);
1480     if (whattodo == 3 || whattodo == 4) tjsch->rctx->reverseonestep = PETSC_TRUE;
1481     if (whattodo == 1) tjsch->rctx->stepsleft = tjsch->rctx->capo-tjsch->rctx->oldcapo;
1482     if (tj->monitor) {
1483       ierr = PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr);
1484       ierr = PetscViewerASCIIPrintf(tj->monitor,"Skip the step from %D to %D (stage values already checkpointed)\n",tjsch->rctx->oldcapo,tjsch->rctx->oldcapo+1);CHKERRQ(ierr);
1485       ierr = PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr);
1486     }
1487     if (!tjsch->rctx->reverseonestep && tjsch->rctx->stepsleft > 0) tjsch->rctx->stepsleft--;
1488     restart++; /* skip one step */
1489   }
1490   if (stack->solution_only || (!stack->solution_only && restart < stepnum)) {
1491     tjsch->recompute = PETSC_TRUE;
1492     ierr = TurnForward(ts);CHKERRQ(ierr);
1493     ierr = ReCompute(ts,tjsch,restart,stepnum);CHKERRQ(ierr);
1494   }
1495   if (!ondisk && ( (stack->solution_only && e->stepnum+1 == stepnum) || (!stack->solution_only && e->stepnum == stepnum) )) {
1496     ierr = StackPop(stack,&e);CHKERRQ(ierr);
1497     ierr = ElementDestroy(stack,e);CHKERRQ(ierr);
1498   }
1499   tjsch->rctx->reverseonestep = PETSC_FALSE;
1500   PetscFunctionReturn(0);
1501 }
1502 #endif
1503 
1504 static PetscErrorCode TSTrajectorySet_Memory(TSTrajectory tj,TS ts,PetscInt stepnum,PetscReal time,Vec X)
1505 {
1506   TJScheduler *tjsch = (TJScheduler*)tj->data;
1507   PetscErrorCode ierr;
1508 
1509   PetscFunctionBegin;
1510   if (!tjsch->recompute) { /* use global stepnum in the forward sweep */
1511     ierr = TSGetTotalSteps(ts,&stepnum);CHKERRQ(ierr);
1512   }
1513   /* for consistency */
1514   if (!tjsch->recompute && stepnum == 0) ts->ptime_prev = ts->ptime-ts->time_step;
1515   switch (tjsch->stype) {
1516     case NONE:
1517       ierr = SetTrajN(ts,tjsch,stepnum,time,X);CHKERRQ(ierr);
1518       break;
1519     case TWO_LEVEL_NOREVOLVE:
1520       ierr = SetTrajTLNR(tj,ts,tjsch,stepnum,time,X);CHKERRQ(ierr);
1521       break;
1522 #if defined(PETSC_HAVE_REVOLVE)
1523     case TWO_LEVEL_REVOLVE:
1524       ierr = SetTrajTLR(tj,ts,tjsch,stepnum,time,X);CHKERRQ(ierr);
1525       break;
1526     case TWO_LEVEL_TWO_REVOLVE:
1527       ierr = SetTrajTLTR(tj,ts,tjsch,stepnum,time,X);CHKERRQ(ierr);
1528       break;
1529     case REVOLVE_OFFLINE:
1530       ierr = SetTrajROF(tj,ts,tjsch,stepnum,time,X);CHKERRQ(ierr);
1531       break;
1532     case REVOLVE_ONLINE:
1533       ierr = SetTrajRON(tj,ts,tjsch,stepnum,time,X);CHKERRQ(ierr);
1534       break;
1535     case REVOLVE_MULTISTAGE:
1536       ierr = SetTrajRMS(tj,ts,tjsch,stepnum,time,X);CHKERRQ(ierr);
1537       break;
1538 #endif
1539     default:
1540       break;
1541   }
1542   PetscFunctionReturn(0);
1543 }
1544 
1545 static PetscErrorCode TSTrajectoryGet_Memory(TSTrajectory tj,TS ts,PetscInt stepnum,PetscReal *t)
1546 {
1547   TJScheduler *tjsch = (TJScheduler*)tj->data;
1548   PetscErrorCode ierr;
1549 
1550   PetscFunctionBegin;
1551   ierr = TSGetTotalSteps(ts,&stepnum);CHKERRQ(ierr);
1552   if (stepnum == 0) PetscFunctionReturn(0);
1553   switch (tjsch->stype) {
1554     case NONE:
1555       ierr = GetTrajN(ts,tjsch,stepnum);CHKERRQ(ierr);
1556       break;
1557     case TWO_LEVEL_NOREVOLVE:
1558       ierr = GetTrajTLNR(tj,ts,tjsch,stepnum);CHKERRQ(ierr);
1559       break;
1560 #if defined(PETSC_HAVE_REVOLVE)
1561     case TWO_LEVEL_REVOLVE:
1562       ierr = GetTrajTLR(tj,ts,tjsch,stepnum);CHKERRQ(ierr);
1563       break;
1564     case TWO_LEVEL_TWO_REVOLVE:
1565       ierr = GetTrajTLTR(tj,ts,tjsch,stepnum);CHKERRQ(ierr);
1566       break;
1567     case REVOLVE_OFFLINE:
1568       ierr = GetTrajROF(tj,ts,tjsch,stepnum);CHKERRQ(ierr);
1569       break;
1570     case REVOLVE_ONLINE:
1571       ierr = GetTrajRON(tj,ts,tjsch,stepnum);CHKERRQ(ierr);
1572       break;
1573     case REVOLVE_MULTISTAGE:
1574       ierr = GetTrajRMS(tj,ts,tjsch,stepnum);CHKERRQ(ierr);
1575       break;
1576 #endif
1577     default:
1578       break;
1579   }
1580   PetscFunctionReturn(0);
1581 }
1582 
1583 PETSC_UNUSED static PetscErrorCode TSTrajectorySetStride_Memory(TSTrajectory tj,PetscInt stride)
1584 {
1585   TJScheduler *tjsch = (TJScheduler*)tj->data;
1586 
1587   PetscFunctionBegin;
1588   tjsch->stride = stride;
1589   PetscFunctionReturn(0);
1590 }
1591 
1592 PETSC_UNUSED static PetscErrorCode TSTrajectorySetMaxCpsRAM_Memory(TSTrajectory tj,PetscInt max_cps_ram)
1593 {
1594   TJScheduler *tjsch = (TJScheduler*)tj->data;
1595 
1596   PetscFunctionBegin;
1597   tjsch->max_cps_ram = max_cps_ram;
1598   PetscFunctionReturn(0);
1599 }
1600 
1601 PETSC_UNUSED static PetscErrorCode TSTrajectorySetMaxCpsDisk_Memory(TSTrajectory tj,PetscInt max_cps_disk)
1602 {
1603   TJScheduler *tjsch = (TJScheduler*)tj->data;
1604 
1605   PetscFunctionBegin;
1606   tjsch->max_cps_disk = max_cps_disk;
1607   PetscFunctionReturn(0);
1608 }
1609 
1610 #if defined(PETSC_HAVE_REVOLVE)
1611 PETSC_UNUSED static PetscErrorCode TSTrajectorySetRevolveOnline(TSTrajectory tj,PetscBool use_online)
1612 {
1613   TJScheduler *tjsch = (TJScheduler*)tj->data;
1614 
1615   PetscFunctionBegin;
1616   tjsch->use_online = use_online;
1617   PetscFunctionReturn(0);
1618 }
1619 #endif
1620 
1621 PETSC_UNUSED static PetscErrorCode TSTrajectorySetSaveStack(TSTrajectory tj,PetscBool save_stack)
1622 {
1623   TJScheduler *tjsch = (TJScheduler*)tj->data;
1624 
1625   PetscFunctionBegin;
1626   tjsch->save_stack = save_stack;
1627   PetscFunctionReturn(0);
1628 }
1629 
1630 PETSC_UNUSED static PetscErrorCode TSTrajectorySetSolutionOnly(TSTrajectory tj,PetscBool solution_only)
1631 {
1632   TJScheduler *tjsch = (TJScheduler*)tj->data;
1633   Stack       *stack = &tjsch->stack;
1634 
1635   PetscFunctionBegin;
1636   stack->solution_only = solution_only;
1637   PetscFunctionReturn(0);
1638 }
1639 
1640 static PetscErrorCode TSTrajectorySetFromOptions_Memory(PetscOptionItems *PetscOptionsObject,TSTrajectory tj)
1641 {
1642   TJScheduler    *tjsch = (TJScheduler*)tj->data;
1643   PetscErrorCode ierr;
1644 
1645   PetscFunctionBegin;
1646   ierr = PetscOptionsHead(PetscOptionsObject,"Memory based TS trajectory options");CHKERRQ(ierr);
1647   {
1648     ierr = PetscOptionsInt("-ts_trajectory_max_cps_ram","Maximum number of checkpoints in RAM","TSTrajectorySetMaxCpsRAM_Memory",tjsch->max_cps_ram,&tjsch->max_cps_ram,NULL);CHKERRQ(ierr);
1649     ierr = PetscOptionsInt("-ts_trajectory_max_cps_disk","Maximum number of checkpoints on disk","TSTrajectorySetMaxCpsDisk_Memory",tjsch->max_cps_disk,&tjsch->max_cps_disk,NULL);CHKERRQ(ierr);
1650     ierr = PetscOptionsInt("-ts_trajectory_stride","Stride to save checkpoints to file","TSTrajectorySetStride_Memory",tjsch->stride,&tjsch->stride,NULL);CHKERRQ(ierr);
1651 #if defined(PETSC_HAVE_REVOLVE)
1652     ierr = PetscOptionsBool("-ts_trajectory_revolve_online","Trick TS trajectory into using online mode of revolve","TSTrajectorySetRevolveOnline",tjsch->use_online,&tjsch->use_online,NULL);CHKERRQ(ierr);
1653 #endif
1654     ierr = PetscOptionsBool("-ts_trajectory_save_stack","Save all stack to disk","TSTrajectorySetSaveStack",tjsch->save_stack,&tjsch->save_stack,NULL);CHKERRQ(ierr);
1655     ierr = PetscOptionsBool("-ts_trajectory_solution_only","Checkpoint solution only","TSTrajectorySetSolutionOnly",tjsch->stack.solution_only,&tjsch->stack.solution_only,NULL);CHKERRQ(ierr);
1656   }
1657   ierr = PetscOptionsTail();CHKERRQ(ierr);
1658   PetscFunctionReturn(0);
1659 }
1660 
1661 static PetscErrorCode TSTrajectorySetUp_Memory(TSTrajectory tj,TS ts)
1662 {
1663   TJScheduler    *tjsch = (TJScheduler*)tj->data;
1664   Stack          *stack = &tjsch->stack;
1665 #if defined(PETSC_HAVE_REVOLVE)
1666   RevolveCTX     *rctx,*rctx2;
1667   DiskStack      *diskstack = &tjsch->diskstack;
1668   PetscInt       diskblocks;
1669 #endif
1670   PetscInt       numY;
1671   PetscBool      flg;
1672   PetscErrorCode ierr;
1673 
1674   PetscFunctionBegin;
1675   PetscStrcmp(((PetscObject)ts->adapt)->type_name,TSADAPTNONE,&flg);
1676   if (flg) tjsch->total_steps = PetscMin(ts->max_steps,(PetscInt)(ceil(ts->max_time/ts->time_step))); /* fixed time step */
1677   if (tjsch->max_cps_ram > 0) stack->stacksize = tjsch->max_cps_ram;
1678 
1679   if (tjsch->stride > 1) { /* two level mode */
1680     if (tjsch->save_stack && tjsch->max_cps_disk > 1 && tjsch->max_cps_disk <= tjsch->max_cps_ram) SETERRQ(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.");
1681     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 */
1682     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 */
1683     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 */
1684   } else { /* single level mode */
1685     if (flg) { /* fixed time step */
1686       if (tjsch->max_cps_ram >= tjsch->total_steps-1 || tjsch->max_cps_ram < 1) tjsch->stype = NONE; /* checkpoint all */
1687       else tjsch->stype = (tjsch->max_cps_disk>1) ? REVOLVE_MULTISTAGE : REVOLVE_OFFLINE;
1688     } else tjsch->stype = NONE; /* checkpoint all for adaptive time step */
1689 #if defined(PETSC_HAVE_REVOLVE)
1690     if (tjsch->use_online) tjsch->stype = REVOLVE_ONLINE; /* trick into online (for testing purpose only) */
1691 #endif
1692   }
1693 
1694   if (tjsch->stype > TWO_LEVEL_NOREVOLVE) {
1695 #ifndef PETSC_HAVE_REVOLVE
1696     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.");
1697 #else
1698     switch (tjsch->stype) {
1699       case TWO_LEVEL_REVOLVE:
1700         revolve_create_offline(tjsch->stride,tjsch->max_cps_ram);
1701         break;
1702       case TWO_LEVEL_TWO_REVOLVE:
1703         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. */
1704         diskstack->stacksize = diskblocks;
1705         revolve_create_offline(tjsch->stride,tjsch->max_cps_ram);
1706         revolve2_create_offline((tjsch->total_steps+tjsch->stride-1)/tjsch->stride,diskblocks);
1707         ierr = PetscCalloc1(1,&rctx2);CHKERRQ(ierr);
1708         rctx2->snaps_in       = diskblocks;
1709         rctx2->reverseonestep = PETSC_FALSE;
1710         rctx2->check          = 0;
1711         rctx2->oldcapo        = 0;
1712         rctx2->capo           = 0;
1713         rctx2->info           = 2;
1714         rctx2->fine           = (tjsch->total_steps+tjsch->stride-1)/tjsch->stride;
1715         tjsch->rctx2          = rctx2;
1716         diskstack->top        = -1;
1717         ierr = PetscMalloc1(diskstack->stacksize*sizeof(PetscInt),&diskstack->container);CHKERRQ(ierr);
1718         break;
1719       case REVOLVE_OFFLINE:
1720         revolve_create_offline(tjsch->total_steps,tjsch->max_cps_ram);
1721         break;
1722       case REVOLVE_ONLINE:
1723         stack->stacksize = tjsch->max_cps_ram;
1724         revolve_create_online(tjsch->max_cps_ram);
1725         break;
1726       case REVOLVE_MULTISTAGE:
1727         revolve_create_multistage(tjsch->total_steps,tjsch->max_cps_ram+tjsch->max_cps_disk,tjsch->max_cps_ram);
1728         break;
1729       default:
1730         break;
1731     }
1732     ierr = PetscCalloc1(1,&rctx);CHKERRQ(ierr);
1733     rctx->snaps_in       = tjsch->max_cps_ram; /* for theta methods snaps_in=2*max_cps_ram */
1734     rctx->reverseonestep = PETSC_FALSE;
1735     rctx->check          = 0;
1736     rctx->oldcapo        = 0;
1737     rctx->capo           = 0;
1738     rctx->info           = 2;
1739     rctx->fine           = (tjsch->stride > 1) ? tjsch->stride : tjsch->total_steps;
1740     tjsch->rctx          = rctx;
1741     if (tjsch->stype == REVOLVE_ONLINE) rctx->fine = -1;
1742 #endif
1743   } else {
1744     if (tjsch->stype == TWO_LEVEL_NOREVOLVE) stack->stacksize = tjsch->stride-1; /* need tjsch->stride-1 at most */
1745     if (tjsch->stype == NONE) {
1746       if (flg) stack->stacksize = stack->solution_only ? tjsch->total_steps : tjsch->total_steps-1; /* fix time step */
1747       else { /* adaptive time step */
1748         if(tjsch->max_cps_ram == -1) stack->stacksize = ts->max_steps; /* if max_cps_ram is not specified, use maximal allowed number of steps for stack size */
1749         tjsch->total_steps = stack->solution_only ? stack->stacksize:stack->stacksize+1; /* will be updated as time integration advances */
1750       }
1751     }
1752   }
1753 
1754   tjsch->recompute = PETSC_FALSE;
1755   ierr = TSGetStages(ts,&numY,PETSC_IGNORE);CHKERRQ(ierr);
1756   ierr = StackCreate(stack,stack->stacksize,numY);CHKERRQ(ierr);
1757   PetscFunctionReturn(0);
1758 }
1759 
1760 static PetscErrorCode TSTrajectoryDestroy_Memory(TSTrajectory tj)
1761 {
1762   TJScheduler    *tjsch = (TJScheduler*)tj->data;
1763   PetscErrorCode ierr;
1764 
1765   PetscFunctionBegin;
1766   if (tjsch->stype > TWO_LEVEL_NOREVOLVE) {
1767 #if defined(PETSC_HAVE_REVOLVE)
1768     revolve_reset();
1769     if (tjsch->stype == TWO_LEVEL_TWO_REVOLVE) {
1770       revolve2_reset();
1771       ierr = PetscFree(tjsch->diskstack.container);CHKERRQ(ierr);
1772     }
1773 #endif
1774   }
1775   ierr = StackDestroy(&tjsch->stack);CHKERRQ(ierr);
1776 #if defined(PETSC_HAVE_REVOLVE)
1777   if (tjsch->stype > TWO_LEVEL_NOREVOLVE) {
1778     ierr = PetscFree(tjsch->rctx);CHKERRQ(ierr);
1779     ierr = PetscFree(tjsch->rctx2);CHKERRQ(ierr);
1780   }
1781 #endif
1782   ierr = PetscFree(tjsch);CHKERRQ(ierr);
1783   PetscFunctionReturn(0);
1784 }
1785 
1786 /*MC
1787       TSTRAJECTORYMEMORY - Stores each solution of the ODE/ADE in memory
1788 
1789   Level: intermediate
1790 
1791 .seealso:  TSTrajectoryCreate(), TS, TSTrajectorySetType()
1792 
1793 M*/
1794 PETSC_EXTERN PetscErrorCode TSTrajectoryCreate_Memory(TSTrajectory tj,TS ts)
1795 {
1796   TJScheduler    *tjsch;
1797   PetscErrorCode ierr;
1798 
1799   PetscFunctionBegin;
1800   tj->ops->set            = TSTrajectorySet_Memory;
1801   tj->ops->get            = TSTrajectoryGet_Memory;
1802   tj->ops->setup          = TSTrajectorySetUp_Memory;
1803   tj->ops->destroy        = TSTrajectoryDestroy_Memory;
1804   tj->ops->setfromoptions = TSTrajectorySetFromOptions_Memory;
1805 
1806   ierr = PetscCalloc1(1,&tjsch);CHKERRQ(ierr);
1807   tjsch->stype        = NONE;
1808   tjsch->max_cps_ram  = -1; /* -1 indicates that it is not set */
1809   tjsch->max_cps_disk = -1; /* -1 indicates that it is not set */
1810   tjsch->stride       = 0; /* if not zero, two-level checkpointing will be used */
1811 #if defined(PETSC_HAVE_REVOLVE)
1812   tjsch->use_online   = PETSC_FALSE;
1813 #endif
1814   tjsch->save_stack   = PETSC_TRUE;
1815 
1816   tjsch->stack.solution_only = PETSC_TRUE;
1817 
1818   tj->data = tjsch;
1819   PetscFunctionReturn(0);
1820 }
1821