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