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