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