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