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