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