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