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