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