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