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