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