xref: /petsc/src/ts/trajectory/impls/memory/trajmemory.c (revision 7409eceac6e206dac4389193de3ec2bf7d26a8ba)
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 
6 /* Limit Revolve to 32-bits */
7 #define PETSC_REVOLVE_INT_MAX  2147483647
8 
9 typedef int PetscRevolveInt;
10 
11 PETSC_STATIC_INLINE PetscErrorCode PetscRevolveIntCast(PetscInt a,PetscRevolveInt *b)
12 {
13   PetscFunctionBegin;
14 #if defined(PETSC_USE_64BIT_INDICES)
15   *b = 0;
16   if (a > PETSC_REVOLVE_INT_MAX) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Parameter is too large for Revolve, which is restricted to 32 bit integers");
17 #endif
18   *b = (PetscRevolveInt)(a);
19   PetscFunctionReturn(0);
20 }
21 #endif
22 #if defined(PETSC_HAVE_CAMS)
23 #include <offline_schedule.h>
24 #endif
25 
26 PetscLogEvent TSTrajectory_DiskWrite, TSTrajectory_DiskRead;
27 static PetscErrorCode TSTrajectorySet_Memory(TSTrajectory,TS,PetscInt,PetscReal,Vec);
28 
29 typedef enum {NONE,TWO_LEVEL_NOREVOLVE,TWO_LEVEL_REVOLVE,TWO_LEVEL_TWO_REVOLVE,REVOLVE_OFFLINE,REVOLVE_ONLINE,REVOLVE_MULTISTAGE,CAMS_OFFLINE} SchedulerType;
30 
31 typedef enum {UNSET=-1,SOLUTIONONLY=0,STAGESONLY=1,SOLUTION_STAGES=2} CheckpointType;
32 
33 #define HaveSolution(m) ((m) == SOLUTIONONLY || (m) == SOLUTION_STAGES)
34 #define HaveStages(m)   ((m) == STAGESONLY || (m) == SOLUTION_STAGES)
35 
36 typedef struct _StackElement {
37   PetscInt       stepnum;
38   Vec            X;
39   Vec            *Y;
40   PetscReal      time;
41   PetscReal      timeprev; /* for no solution_only mode */
42   PetscReal      timenext; /* for solution_only mode */
43   CheckpointType cptype;
44 } *StackElement;
45 
46 #if defined(PETSC_HAVE_REVOLVE)
47 typedef struct _RevolveCTX {
48   PetscBool       reverseonestep;
49   PetscRevolveInt where;
50   PetscRevolveInt snaps_in;
51   PetscRevolveInt stepsleft;
52   PetscRevolveInt check;
53   PetscRevolveInt oldcapo;
54   PetscRevolveInt capo;
55   PetscRevolveInt fine;
56   PetscRevolveInt info;
57 } RevolveCTX;
58 #endif
59 
60 #if defined(PETSC_HAVE_CAMS)
61 typedef struct _CAMSCTX {
62   PetscInt lastcheckpointstep;
63   PetscInt lastcheckpointtype;
64   PetscInt num_units_avail;
65   PetscInt endstep;
66   PetscInt num_stages;
67   PetscInt nextcheckpointstep;
68   PetscInt nextcheckpointtype; /* (0) solution only (1) stages (2) solution+stages */
69   PetscInt info;
70 } CAMSCTX;
71 #endif
72 
73 typedef struct _Stack {
74   PetscInt      stacksize;
75   PetscInt      top;
76   StackElement  *container;
77   PetscInt      nallocated;
78   PetscInt      numY;
79   PetscBool     solution_only;
80   PetscBool     use_dram;
81 } Stack;
82 
83 typedef struct _DiskStack {
84   PetscInt  stacksize;
85   PetscInt  top;
86   PetscInt  *container;
87 } DiskStack;
88 
89 typedef struct _TJScheduler {
90   SchedulerType stype;
91 #if defined(PETSC_HAVE_REVOLVE)
92   RevolveCTX    *rctx,*rctx2;
93   PetscBool     use_online;
94   PetscInt      store_stride;
95 #endif
96 #if defined(PETSC_HAVE_CAMS)
97   CAMSCTX       *actx;
98 #endif
99   PetscBool     recompute;
100   PetscBool     skip_trajectory;
101   PetscBool     save_stack;
102   PetscInt      max_units_ram;  /* maximum checkpointing units in RAM */
103   PetscInt      max_units_disk; /* maximum checkpointing units on disk */
104   PetscInt      max_cps_ram;    /* maximum checkpoints in RAM */
105   PetscInt      max_cps_disk;   /* maximum checkpoints on disk */
106   PetscInt      stride;
107   PetscInt      total_steps;    /* total number of steps */
108   Stack         stack;
109   DiskStack     diskstack;
110   PetscViewer   viewer;
111 } TJScheduler;
112 
113 static PetscErrorCode TurnForwardWithStepsize(TS ts,PetscReal nextstepsize)
114 {
115   PetscErrorCode ierr;
116 
117   PetscFunctionBegin;
118   /* reverse the direction */
119   ierr = TSSetTimeStep(ts,nextstepsize);CHKERRQ(ierr);
120   PetscFunctionReturn(0);
121 }
122 
123 static PetscErrorCode TurnForward(TS ts)
124 {
125   PetscReal      stepsize;
126   PetscErrorCode ierr;
127 
128   PetscFunctionBegin;
129   /* reverse the direction */
130   ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr);
131   ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr);
132   PetscFunctionReturn(0);
133 }
134 
135 static PetscErrorCode TurnBackward(TS ts)
136 {
137   PetscReal      stepsize;
138   PetscErrorCode ierr;
139 
140   PetscFunctionBegin;
141   if (!ts->trajectory->adjoint_solve_mode) PetscFunctionReturn(0);
142   /* reverse the direction */
143   stepsize = ts->ptime_prev-ts->ptime;
144   ierr = TSSetTimeStep(ts,stepsize);CHKERRQ(ierr);
145   PetscFunctionReturn(0);
146 }
147 
148 static PetscErrorCode ElementCreate(TS ts,CheckpointType cptype,Stack *stack,StackElement *e)
149 {
150   Vec            X;
151   Vec            *Y;
152   PetscErrorCode ierr;
153 
154   PetscFunctionBegin;
155   if (stack->top < stack->stacksize-1 && stack->container[stack->top+1]) {
156     *e = stack->container[stack->top+1];
157     if (HaveSolution(cptype) && !(*e)->X) {
158       ierr = TSGetSolution(ts,&X);CHKERRQ(ierr);
159       ierr = VecDuplicate(X,&(*e)->X);CHKERRQ(ierr);
160     }
161     if (cptype==1 && (*e)->X) {
162       ierr = VecDestroy(&(*e)->X);CHKERRQ(ierr);
163     }
164     if (HaveStages(cptype) && !(*e)->Y) {
165       ierr = TSGetStages(ts,&stack->numY,&Y);CHKERRQ(ierr);
166       if (stack->numY) {
167         ierr = VecDuplicateVecs(Y[0],stack->numY,&(*e)->Y);CHKERRQ(ierr);
168       }
169     }
170     if (cptype==0 && (*e)->Y) {
171       ierr = VecDestroyVecs(stack->numY,&(*e)->Y);CHKERRQ(ierr);
172     }
173     (*e)->cptype = cptype;
174     PetscFunctionReturn(0);
175   }
176   if (stack->use_dram) {
177     ierr = PetscMallocSetDRAM();CHKERRQ(ierr);
178   }
179   ierr = PetscNew(e);CHKERRQ(ierr);
180   if (HaveSolution(cptype)) {
181     ierr = TSGetSolution(ts,&X);CHKERRQ(ierr);
182     ierr = VecDuplicate(X,&(*e)->X);CHKERRQ(ierr);
183   }
184   if (HaveStages(cptype)) {
185     ierr = TSGetStages(ts,&stack->numY,&Y);CHKERRQ(ierr);
186     if (stack->numY) {
187       ierr = VecDuplicateVecs(Y[0],stack->numY,&(*e)->Y);CHKERRQ(ierr);
188     }
189   }
190   if (stack->use_dram) {
191     ierr = PetscMallocResetDRAM();CHKERRQ(ierr);
192   }
193   stack->nallocated++;
194   (*e)->cptype = cptype;
195   PetscFunctionReturn(0);
196 }
197 
198 static PetscErrorCode ElementSet(TS ts, Stack *stack, StackElement *e, PetscInt stepnum, PetscReal time, Vec X)
199 {
200   Vec            *Y;
201   PetscInt       i;
202   PetscReal      timeprev;
203   PetscErrorCode ierr;
204 
205   PetscFunctionBegin;
206   if (HaveSolution((*e)->cptype)) {
207     ierr = VecCopy(X,(*e)->X);CHKERRQ(ierr);
208   }
209   if (HaveStages((*e)->cptype)) {
210     ierr = TSGetStages(ts,&stack->numY,&Y);CHKERRQ(ierr);
211     for (i=0;i<stack->numY;i++) {
212       ierr = VecCopy(Y[i],(*e)->Y[i]);CHKERRQ(ierr);
213     }
214   }
215   (*e)->stepnum = stepnum;
216   (*e)->time    = time;
217   /* for consistency */
218   if (stepnum == 0) {
219     (*e)->timeprev = (*e)->time - ts->time_step;
220   } else {
221     ierr = TSGetPrevTime(ts,&timeprev);CHKERRQ(ierr);
222     (*e)->timeprev = timeprev;
223   }
224   PetscFunctionReturn(0);
225 }
226 
227 static PetscErrorCode ElementDestroy(Stack *stack,StackElement e)
228 {
229   PetscErrorCode ierr;
230 
231   PetscFunctionBegin;
232   if (stack->use_dram) {
233     ierr = PetscMallocSetDRAM();CHKERRQ(ierr);
234   }
235   ierr = VecDestroy(&e->X);CHKERRQ(ierr);
236   if (e->Y) {
237     ierr = VecDestroyVecs(stack->numY,&e->Y);CHKERRQ(ierr);
238   }
239   ierr = PetscFree(e);CHKERRQ(ierr);
240   if (stack->use_dram) {
241     ierr = PetscMallocResetDRAM();CHKERRQ(ierr);
242   }
243   stack->nallocated--;
244   PetscFunctionReturn(0);
245 }
246 
247 static PetscErrorCode StackResize(Stack *stack,PetscInt newsize)
248 {
249   StackElement   *newcontainer;
250   PetscInt       i;
251   PetscErrorCode ierr;
252 
253   PetscFunctionBegin;
254   ierr = PetscCalloc1(newsize*sizeof(StackElement),&newcontainer);CHKERRQ(ierr);
255   for (i=0;i<stack->stacksize;i++) {
256     newcontainer[i] = stack->container[i];
257   }
258   ierr = PetscFree(stack->container);CHKERRQ(ierr);
259   stack->container = newcontainer;
260   stack->stacksize = newsize;
261   PetscFunctionReturn(0);
262 }
263 
264 static PetscErrorCode StackPush(Stack *stack,StackElement e)
265 {
266   PetscFunctionBegin;
267   if (stack->top+1 >= stack->stacksize) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MEMC,"Maximum stack size (%D) exceeded",stack->stacksize);
268   stack->container[++stack->top] = e;
269   PetscFunctionReturn(0);
270 }
271 
272 static PetscErrorCode StackPop(Stack *stack,StackElement *e)
273 {
274   PetscFunctionBegin;
275   *e = NULL;
276   if (stack->top == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEMC,"Empty stack");
277   *e = stack->container[stack->top--];
278   PetscFunctionReturn(0);
279 }
280 
281 static PetscErrorCode StackTop(Stack *stack,StackElement *e)
282 {
283   PetscFunctionBegin;
284   *e = stack->container[stack->top];
285   PetscFunctionReturn(0);
286 }
287 
288 static PetscErrorCode StackInit(Stack *stack,PetscInt size,PetscInt ny)
289 {
290   PetscErrorCode ierr;
291 
292   PetscFunctionBegin;
293   stack->top  = -1;
294   stack->numY = ny;
295 
296   if (!stack->container) {
297     ierr = PetscCalloc1(size,&stack->container);CHKERRQ(ierr);
298   }
299   PetscFunctionReturn(0);
300 }
301 
302 static PetscErrorCode StackDestroy(Stack *stack)
303 {
304   PetscInt       i,n = stack->nallocated;
305   PetscErrorCode ierr;
306 
307   PetscFunctionBegin;
308   if (!stack->container) PetscFunctionReturn(0);
309   if (stack->top+1 > stack->nallocated) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Stack size does not match element counter %D",stack->nallocated);
310   for (i=0; i<n; i++) {
311     ierr = ElementDestroy(stack,stack->container[i]);CHKERRQ(ierr);
312   }
313   ierr = PetscFree(stack->container);CHKERRQ(ierr);
314   PetscFunctionReturn(0);
315 }
316 
317 static PetscErrorCode StackFind(Stack *stack,StackElement *e,PetscInt index)
318 {
319   PetscFunctionBegin;
320   *e = NULL;
321   if (index < 0 || index > stack->top) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Invalid index %D",index);
322   *e = stack->container[index];
323   PetscFunctionReturn(0);
324 }
325 
326 static PetscErrorCode WriteToDisk(PetscInt stepnum,PetscReal time,PetscReal timeprev,Vec X,Vec *Y,PetscInt numY,CheckpointType cptype,PetscViewer viewer)
327 {
328   PetscInt       i;
329   PetscErrorCode ierr;
330 
331   PetscFunctionBegin;
332   ierr = PetscViewerBinaryWrite(viewer,&stepnum,1,PETSC_INT);CHKERRQ(ierr);
333   if (HaveSolution(cptype)) {
334     ierr = VecView(X,viewer);CHKERRQ(ierr);
335   }
336   if (HaveStages(cptype)) {
337     for (i=0;i<numY;i++) {
338       ierr = VecView(Y[i],viewer);CHKERRQ(ierr);
339     }
340   }
341   ierr = PetscViewerBinaryWrite(viewer,&time,1,PETSC_REAL);CHKERRQ(ierr);
342   ierr = PetscViewerBinaryWrite(viewer,&timeprev,1,PETSC_REAL);CHKERRQ(ierr);
343   PetscFunctionReturn(0);
344 }
345 
346 static PetscErrorCode ReadFromDisk(PetscInt *stepnum,PetscReal *time,PetscReal *timeprev,Vec X,Vec *Y,PetscInt numY,CheckpointType cptype,PetscViewer viewer)
347 {
348   PetscInt       i;
349   PetscErrorCode ierr;
350 
351   PetscFunctionBegin;
352   ierr = PetscViewerBinaryRead(viewer,stepnum,1,NULL,PETSC_INT);CHKERRQ(ierr);
353   if (HaveSolution(cptype)) {
354     ierr = VecLoad(X,viewer);CHKERRQ(ierr);
355   }
356   if (HaveStages(cptype)) {
357     for (i=0;i<numY;i++) {
358       ierr = VecLoad(Y[i],viewer);CHKERRQ(ierr);
359     }
360   }
361   ierr = PetscViewerBinaryRead(viewer,time,1,NULL,PETSC_REAL);CHKERRQ(ierr);
362   ierr = PetscViewerBinaryRead(viewer,timeprev,1,NULL,PETSC_REAL);CHKERRQ(ierr);
363   PetscFunctionReturn(0);
364 }
365 
366 static PetscErrorCode StackDumpAll(TSTrajectory tj,TS ts,Stack *stack,PetscInt id)
367 {
368   Vec            *Y;
369   PetscInt       i,ndumped,cptype_int;
370   StackElement   e = NULL;
371   TJScheduler    *tjsch = (TJScheduler*)tj->data;
372   char           filename[PETSC_MAX_PATH_LEN];
373   PetscErrorCode ierr;
374   MPI_Comm       comm;
375 
376   PetscFunctionBegin;
377   ierr = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr);
378   if (tj->monitor) {
379     ierr = PetscViewerASCIIPushTab(tj->monitor);CHKERRQ(ierr);
380     ierr = PetscViewerASCIIPrintf(tj->monitor,"Dump stack id %D to file\n",id);CHKERRQ(ierr);
381     ierr = PetscViewerASCIIPopTab(tj->monitor);CHKERRQ(ierr);
382   }
383   ierr = PetscSNPrintf(filename,sizeof(filename),"%s/TS-STACK%06d.bin",tj->dirname,id);CHKERRQ(ierr);
384   ierr = PetscViewerFileSetName(tjsch->viewer,filename);CHKERRQ(ierr);
385   ierr = PetscViewerSetUp(tjsch->viewer);CHKERRQ(ierr);
386   ndumped = stack->top+1;
387   ierr = PetscViewerBinaryWrite(tjsch->viewer,&ndumped,1,PETSC_INT);CHKERRQ(ierr);
388   for (i=0;i<ndumped;i++) {
389     e = stack->container[i];
390     cptype_int = (PetscInt)e->cptype;
391     ierr = PetscViewerBinaryWrite(tjsch->viewer,&cptype_int,1,PETSC_INT);CHKERRQ(ierr);
392     ierr = PetscLogEventBegin(TSTrajectory_DiskWrite,tj,ts,0,0);CHKERRQ(ierr);
393     ierr = WriteToDisk(e->stepnum,e->time,e->timeprev,e->X,e->Y,stack->numY,e->cptype,tjsch->viewer);CHKERRQ(ierr);
394     ierr = PetscLogEventEnd(TSTrajectory_DiskWrite,tj,ts,0,0);CHKERRQ(ierr);
395     ts->trajectory->diskwrites++;
396     ierr = StackPop(stack,&e);CHKERRQ(ierr);
397   }
398   /* 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 */
399   ierr = TSGetStages(ts,&stack->numY,&Y);CHKERRQ(ierr);
400   ierr = PetscLogEventBegin(TSTrajectory_DiskWrite,tj,ts,0,0);CHKERRQ(ierr);
401   ierr = WriteToDisk(ts->steps,ts->ptime,ts->ptime_prev,ts->vec_sol,Y,stack->numY,SOLUTION_STAGES,tjsch->viewer);CHKERRQ(ierr);
402   ierr = PetscLogEventEnd(TSTrajectory_DiskWrite,tj,ts,0,0);CHKERRQ(ierr);
403   ts->trajectory->diskwrites++;
404   PetscFunctionReturn(0);
405 }
406 
407 static PetscErrorCode StackLoadAll(TSTrajectory tj,TS ts,Stack *stack,PetscInt id)
408 {
409   Vec            *Y;
410   PetscInt       i,nloaded,cptype_int;
411   StackElement   e;
412   PetscViewer    viewer;
413   char           filename[PETSC_MAX_PATH_LEN];
414   PetscErrorCode ierr;
415 
416   PetscFunctionBegin;
417   if (tj->monitor) {
418     ierr = PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr);
419     ierr = PetscViewerASCIIPrintf(tj->monitor,"Load stack from file\n");CHKERRQ(ierr);
420     ierr = PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr);
421   }
422   ierr = PetscSNPrintf(filename,sizeof filename,"%s/TS-STACK%06d.bin",tj->dirname,id);CHKERRQ(ierr);
423   ierr = PetscViewerBinaryOpen(PetscObjectComm((PetscObject)tj),filename,FILE_MODE_READ,&viewer);CHKERRQ(ierr);
424   ierr = PetscViewerBinarySetSkipInfo(viewer,PETSC_TRUE);CHKERRQ(ierr);
425   ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_NATIVE);CHKERRQ(ierr);
426   ierr = PetscViewerBinaryRead(viewer,&nloaded,1,NULL,PETSC_INT);CHKERRQ(ierr);
427   for (i=0;i<nloaded;i++) {
428     ierr = PetscViewerBinaryRead(viewer,&cptype_int,1,NULL,PETSC_INT);CHKERRQ(ierr);
429     ierr = ElementCreate(ts,(CheckpointType)cptype_int,stack,&e);CHKERRQ(ierr);
430     ierr = StackPush(stack,e);CHKERRQ(ierr);
431     ierr = PetscLogEventBegin(TSTrajectory_DiskRead,tj,ts,0,0);CHKERRQ(ierr);
432     ierr = ReadFromDisk(&e->stepnum,&e->time,&e->timeprev,e->X,e->Y,stack->numY,e->cptype,viewer);CHKERRQ(ierr);
433     ierr = PetscLogEventEnd(TSTrajectory_DiskRead,tj,ts,0,0);CHKERRQ(ierr);
434     ts->trajectory->diskreads++;
435   }
436   /* load the last step into TS */
437   ierr = TSGetStages(ts,&stack->numY,&Y);CHKERRQ(ierr);
438   ierr = PetscLogEventBegin(TSTrajectory_DiskRead,tj,ts,0,0);CHKERRQ(ierr);
439   ierr = ReadFromDisk(&ts->steps,&ts->ptime,&ts->ptime_prev,ts->vec_sol,Y,stack->numY,SOLUTION_STAGES,viewer);CHKERRQ(ierr);
440   ierr = PetscLogEventEnd(TSTrajectory_DiskRead,tj,ts,0,0);CHKERRQ(ierr);
441   ts->trajectory->diskreads++;
442   ierr = TurnBackward(ts);CHKERRQ(ierr);
443   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
444   PetscFunctionReturn(0);
445 }
446 
447 #if defined(PETSC_HAVE_REVOLVE)
448 static PetscErrorCode StackLoadLast(TSTrajectory tj,TS ts,Stack *stack,PetscInt id)
449 {
450   Vec            *Y;
451   PetscInt       size;
452   PetscViewer    viewer;
453   char           filename[PETSC_MAX_PATH_LEN];
454 #if defined(PETSC_HAVE_MPIIO)
455   PetscBool      usempiio;
456 #endif
457   int            fd;
458   off_t          off,offset;
459   PetscErrorCode ierr;
460 
461   PetscFunctionBegin;
462   if (tj->monitor) {
463     ierr = PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr);
464     ierr = PetscViewerASCIIPrintf(tj->monitor,"Load last stack element from file\n");CHKERRQ(ierr);
465     ierr = PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr);
466   }
467   ierr = TSGetStages(ts,&stack->numY,&Y);CHKERRQ(ierr);
468   ierr = VecGetSize(Y[0],&size);CHKERRQ(ierr);
469   /* VecView writes to file two extra int's for class id and number of rows */
470   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;
471 
472   ierr = PetscSNPrintf(filename,sizeof filename,"%s/TS-STACK%06d.bin",tj->dirname,id);CHKERRQ(ierr);
473   ierr = PetscViewerBinaryOpen(PetscObjectComm((PetscObject)tj),filename,FILE_MODE_READ,&viewer);CHKERRQ(ierr);
474   ierr = PetscViewerBinarySetSkipInfo(viewer,PETSC_TRUE);CHKERRQ(ierr);
475   ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_NATIVE);CHKERRQ(ierr);
476 #if defined(PETSC_HAVE_MPIIO)
477   ierr = PetscViewerBinaryGetUseMPIIO(viewer,&usempiio);CHKERRQ(ierr);
478   if (usempiio) {
479     ierr = PetscViewerBinaryGetMPIIODescriptor(viewer,(MPI_File*)&fd);CHKERRQ(ierr);
480     ierr = PetscBinarySynchronizedSeek(PetscObjectComm((PetscObject)tj),fd,off,PETSC_BINARY_SEEK_END,&offset);CHKERRQ(ierr);
481   } else {
482 #endif
483     ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
484     ierr = PetscBinarySeek(fd,off,PETSC_BINARY_SEEK_END,&offset);CHKERRQ(ierr);
485 #if defined(PETSC_HAVE_MPIIO)
486   }
487 #endif
488   /* load the last step into TS */
489   ierr = PetscLogEventBegin(TSTrajectory_DiskRead,tj,ts,0,0);CHKERRQ(ierr);
490   ierr = ReadFromDisk(&ts->steps,&ts->ptime,&ts->ptime_prev,ts->vec_sol,Y,stack->numY,SOLUTION_STAGES,viewer);CHKERRQ(ierr);
491   ierr = PetscLogEventEnd(TSTrajectory_DiskRead,tj,ts,0,0);CHKERRQ(ierr);
492   ts->trajectory->diskreads++;
493   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
494   ierr = TurnBackward(ts);CHKERRQ(ierr);
495   PetscFunctionReturn(0);
496 }
497 #endif
498 
499 static PetscErrorCode DumpSingle(TSTrajectory tj,TS ts,Stack *stack,PetscInt id)
500 {
501   Vec            *Y;
502   PetscInt       stepnum;
503   TJScheduler    *tjsch = (TJScheduler*)tj->data;
504   char           filename[PETSC_MAX_PATH_LEN];
505   PetscErrorCode ierr;
506   MPI_Comm       comm;
507 
508   PetscFunctionBegin;
509   ierr = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr);
510   if (tj->monitor) {
511     ierr = PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr);
512     ierr = PetscViewerASCIIPrintf(tj->monitor,"Dump a single point from file\n");CHKERRQ(ierr);
513     ierr = PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr);
514   }
515   ierr = TSGetStepNumber(ts,&stepnum);CHKERRQ(ierr);
516   ierr = PetscSNPrintf(filename,sizeof(filename),"%s/TS-CPS%06d.bin",tj->dirname,id);CHKERRQ(ierr);
517   ierr = PetscViewerFileSetName(tjsch->viewer,filename);CHKERRQ(ierr);
518   ierr = PetscViewerSetUp(tjsch->viewer);CHKERRQ(ierr);
519 
520   ierr = TSGetStages(ts,&stack->numY,&Y);CHKERRQ(ierr);
521   ierr = PetscLogEventBegin(TSTrajectory_DiskWrite,tj,ts,0,0);CHKERRQ(ierr);
522   ierr = WriteToDisk(stepnum,ts->ptime,ts->ptime_prev,ts->vec_sol,Y,stack->numY,SOLUTION_STAGES,tjsch->viewer);CHKERRQ(ierr);
523   ierr = PetscLogEventEnd(TSTrajectory_DiskWrite,tj,ts,0,0);CHKERRQ(ierr);
524   ts->trajectory->diskwrites++;
525   PetscFunctionReturn(0);
526 }
527 
528 static PetscErrorCode LoadSingle(TSTrajectory tj,TS ts,Stack *stack,PetscInt id)
529 {
530   Vec            *Y;
531   PetscViewer    viewer;
532   char           filename[PETSC_MAX_PATH_LEN];
533   PetscErrorCode ierr;
534 
535   PetscFunctionBegin;
536   if (tj->monitor) {
537     ierr = PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr);
538     ierr = PetscViewerASCIIPrintf(tj->monitor,"Load a single point from file\n");CHKERRQ(ierr);
539     ierr = PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr);
540   }
541   ierr = PetscSNPrintf(filename,sizeof filename,"%s/TS-CPS%06d.bin",tj->dirname,id);CHKERRQ(ierr);
542   ierr = PetscViewerBinaryOpen(PetscObjectComm((PetscObject)tj),filename,FILE_MODE_READ,&viewer);CHKERRQ(ierr);
543   ierr = PetscViewerBinarySetSkipInfo(viewer,PETSC_TRUE);CHKERRQ(ierr);
544   ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_NATIVE);CHKERRQ(ierr);
545   ierr = TSGetStages(ts,&stack->numY,&Y);CHKERRQ(ierr);
546   ierr = PetscLogEventBegin(TSTrajectory_DiskRead,tj,ts,0,0);CHKERRQ(ierr);
547   ierr = ReadFromDisk(&ts->steps,&ts->ptime,&ts->ptime_prev,ts->vec_sol,Y,stack->numY,SOLUTION_STAGES,viewer);CHKERRQ(ierr);
548   ierr = PetscLogEventEnd(TSTrajectory_DiskRead,tj,ts,0,0);CHKERRQ(ierr);
549   ts->trajectory->diskreads++;
550   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
551   PetscFunctionReturn(0);
552 }
553 
554 static PetscErrorCode UpdateTS(TS ts,Stack *stack,StackElement e,PetscInt stepnum,PetscBool adjoint_mode)
555 {
556   Vec            *Y;
557   PetscInt       i;
558   PetscErrorCode ierr;
559 
560   PetscFunctionBegin;
561   if (HaveSolution(e->cptype) && e->stepnum!=stepnum) {
562     ierr = VecCopy(e->X,ts->vec_sol);CHKERRQ(ierr);
563   }
564   if (HaveStages(e->cptype)) {
565     ierr = TSGetStages(ts,&stack->numY,&Y);CHKERRQ(ierr);
566     if (e->stepnum && e->stepnum==stepnum) {
567       for (i=0;i<stack->numY;i++) {
568         ierr = VecCopy(e->Y[i],Y[i]);CHKERRQ(ierr);
569       }
570     } else if (ts->stifflyaccurate) {
571       ierr = VecCopy(e->Y[stack->numY-1],ts->vec_sol);CHKERRQ(ierr);
572     }
573   }
574   if (adjoint_mode) {
575     ierr = TSSetTimeStep(ts,e->timeprev-e->time);CHKERRQ(ierr); /* stepsize will be negative */
576   } else {
577     ierr = TSSetTimeStep(ts,e->time-e->timeprev);CHKERRQ(ierr); /* stepsize will be positive */
578   }
579   ts->ptime      = e->time;
580   ts->ptime_prev = e->timeprev;
581   PetscFunctionReturn(0);
582 }
583 
584 static PetscErrorCode ReCompute(TS ts,TJScheduler *tjsch,PetscInt stepnumbegin,PetscInt stepnumend)
585 {
586   Stack          *stack = &tjsch->stack;
587   PetscInt       i;
588   PetscErrorCode ierr;
589 
590   PetscFunctionBegin;
591   tjsch->recompute = PETSC_TRUE; /* hints TSTrajectorySet() that it is in recompute mode */
592   ierr = TSSetStepNumber(ts,stepnumbegin);CHKERRQ(ierr);/* global step number */
593   for (i=stepnumbegin;i<stepnumend;i++) { /* assume fixed step size */
594     if (stack->solution_only && !tjsch->skip_trajectory) { /* revolve online need this */
595       /* don't use the public interface as it will update the TSHistory: this need a better fix */
596       ierr = TSTrajectorySet_Memory(ts->trajectory,ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
597     }
598     ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
599     ierr = TSStep(ts);CHKERRQ(ierr);
600     if (!stack->solution_only && !tjsch->skip_trajectory) {
601       /* don't use the public interface as it will update the TSHistory: this need a better fix */
602       ierr = TSTrajectorySet_Memory(ts->trajectory,ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
603     }
604     ierr = TSEventHandler(ts);CHKERRQ(ierr);
605     if (!ts->steprollback) {
606       ierr = TSPostStep(ts);CHKERRQ(ierr);
607     }
608   }
609   ierr = TurnBackward(ts);CHKERRQ(ierr);
610   ts->trajectory->recomps += stepnumend-stepnumbegin; /* recomputation counter */
611   ierr = TSSetStepNumber(ts,stepnumend);CHKERRQ(ierr);
612   tjsch->recompute = PETSC_FALSE; /* reset the flag for recompute mode */
613   PetscFunctionReturn(0);
614 }
615 
616 static PetscErrorCode TopLevelStore(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscInt localstepnum,PetscInt laststridesize,PetscBool *done)
617 {
618   Stack          *stack = &tjsch->stack;
619   DiskStack      *diskstack = &tjsch->diskstack;
620   PetscInt       stridenum;
621   PetscErrorCode ierr;
622 
623   PetscFunctionBegin;
624   *done = PETSC_FALSE;
625   stridenum    = stepnum/tjsch->stride;
626   /* make sure saved checkpoint id starts from 1
627      skip last stride when using stridenum+1
628      skip first stride when using stridenum */
629   if (stack->solution_only) {
630     if (tjsch->save_stack) {
631       if (localstepnum == tjsch->stride-1 && stepnum < tjsch->total_steps-laststridesize) { /* current step will be saved without going through stack */
632         ierr = StackDumpAll(tj,ts,stack,stridenum+1);CHKERRQ(ierr);
633         if (tjsch->stype == TWO_LEVEL_TWO_REVOLVE) diskstack->container[++diskstack->top] = stridenum+1;
634         *done = PETSC_TRUE;
635       }
636     } else {
637       if (localstepnum == 0 && stepnum < tjsch->total_steps-laststridesize) {
638         ierr = DumpSingle(tj,ts,stack,stridenum+1);CHKERRQ(ierr);
639         if (tjsch->stype == TWO_LEVEL_TWO_REVOLVE) diskstack->container[++diskstack->top] = stridenum+1;
640         *done = PETSC_TRUE;
641       }
642     }
643   } else {
644     if (tjsch->save_stack) {
645       if (localstepnum == 0 && stepnum < tjsch->total_steps && stepnum != 0) { /* skip the first stride */
646         ierr = StackDumpAll(tj,ts,stack,stridenum);CHKERRQ(ierr);
647         if (tjsch->stype == TWO_LEVEL_TWO_REVOLVE) diskstack->container[++diskstack->top] = stridenum;
648         *done = PETSC_TRUE;
649       }
650     } else {
651       if (localstepnum == 1 && stepnum < tjsch->total_steps-laststridesize) {
652         ierr = DumpSingle(tj,ts,stack,stridenum+1);CHKERRQ(ierr);
653         if (tjsch->stype == TWO_LEVEL_TWO_REVOLVE) diskstack->container[++diskstack->top] = stridenum+1;
654         *done = PETSC_TRUE;
655       }
656     }
657   }
658   PetscFunctionReturn(0);
659 }
660 
661 static PetscErrorCode SetTrajN(TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X)
662 {
663   Stack          *stack = &tjsch->stack;
664   StackElement   e;
665   CheckpointType cptype;
666   PetscErrorCode ierr;
667 
668   PetscFunctionBegin;
669   /* skip the last step */
670   if (ts->reason) { /* only affect the forward run */
671     /* update total_steps in the end of forward run */
672     if (stepnum != tjsch->total_steps) tjsch->total_steps = stepnum;
673     if (stack->solution_only) {
674       /* get rid of the solution at second last step */
675       ierr = StackPop(stack,&e);CHKERRQ(ierr);
676     }
677     PetscFunctionReturn(0);
678   }
679   /*  do not save trajectory at the recompute stage for solution_only mode */
680   if (stack->solution_only && tjsch->recompute) PetscFunctionReturn(0);
681   /* skip the first step for no_solution_only mode */
682   if (!stack->solution_only && stepnum == 0) PetscFunctionReturn(0);
683 
684   /* resize the stack */
685   if (stack->top+1 == stack->stacksize) {
686     ierr = StackResize(stack,2*stack->stacksize);CHKERRQ(ierr);
687   }
688   /* update timenext for the previous step; necessary for step adaptivity */
689   if (stack->top > -1) {
690     ierr = StackTop(stack,&e);CHKERRQ(ierr);
691     e->timenext = ts->ptime;
692   }
693   if (stepnum < stack->top) {
694     SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_MEMC,"Illegal modification of a non-top stack element");
695   }
696   cptype = stack->solution_only ? SOLUTIONONLY : STAGESONLY;
697   ierr = ElementCreate(ts,cptype,stack,&e);CHKERRQ(ierr);
698   ierr = ElementSet(ts,stack,&e,stepnum,time,X);CHKERRQ(ierr);
699   ierr = StackPush(stack,e);CHKERRQ(ierr);
700   PetscFunctionReturn(0);
701 }
702 
703 static PetscErrorCode SetTrajN_2(TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X)
704 {
705   Stack          *stack = &tjsch->stack;
706   StackElement   e;
707   CheckpointType cptype;
708   PetscErrorCode ierr;
709 
710   PetscFunctionBegin;
711   if (stack->top+1 == stack->stacksize) {
712     ierr = StackResize(stack,2*stack->stacksize);CHKERRQ(ierr);
713   }
714   /* update timenext for the previous step; necessary for step adaptivity */
715   if (stack->top > -1) {
716     ierr = StackTop(stack,&e);CHKERRQ(ierr);
717     e->timenext = ts->ptime;
718   }
719   if (stepnum < stack->top) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_MEMC,"Illegal modification of a non-top stack element");
720   cptype = stack->solution_only ? SOLUTIONONLY : SOLUTION_STAGES; /* Always include solution in a checkpoint in non-adjoint mode */
721   ierr = ElementCreate(ts,cptype,stack,&e);CHKERRQ(ierr);
722   ierr = ElementSet(ts,stack,&e,stepnum,time,X);CHKERRQ(ierr);
723   ierr = StackPush(stack,e);CHKERRQ(ierr);
724   PetscFunctionReturn(0);
725 }
726 
727 static PetscErrorCode GetTrajN(TS ts,TJScheduler *tjsch,PetscInt stepnum)
728 {
729   Stack          *stack = &tjsch->stack;
730   StackElement   e;
731   PetscInt       ns;
732   PetscErrorCode ierr;
733 
734   PetscFunctionBegin;
735   /* If TSTrajectoryGet() is called after TSAdjointSolve() converges (e.g. outside the while loop in TSAdjointSolve()), skip getting the checkpoint. */
736   if (ts->reason) PetscFunctionReturn(0);
737   if (stepnum == tjsch->total_steps) {
738     ierr = TurnBackward(ts);CHKERRQ(ierr);
739     PetscFunctionReturn(0);
740   }
741   /* restore a checkpoint */
742   ierr = StackTop(stack,&e);CHKERRQ(ierr);
743   ierr = UpdateTS(ts,stack,e,stepnum,PETSC_TRUE);CHKERRQ(ierr);
744   ierr = TSGetStages(ts,&ns,PETSC_IGNORE);CHKERRQ(ierr);
745   if (stack->solution_only && ns) { /* recompute one step */
746     ierr = TurnForwardWithStepsize(ts,e->timenext-e->time);CHKERRQ(ierr);
747     ierr = ReCompute(ts,tjsch,e->stepnum,stepnum);CHKERRQ(ierr);
748   }
749   ierr = StackPop(stack,&e);CHKERRQ(ierr);
750   PetscFunctionReturn(0);
751 }
752 
753 static PetscErrorCode GetTrajN_2(TS ts,TJScheduler *tjsch,PetscInt stepnum)
754 {
755   Stack          *stack = &tjsch->stack;
756   StackElement   e = NULL;
757   PetscErrorCode ierr;
758 
759   PetscFunctionBegin;
760   ierr = StackFind(stack,&e,stepnum);CHKERRQ(ierr);
761   if (stepnum != e->stepnum) SETERRQ2(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Inconsistent steps! %D != %D",stepnum,e->stepnum);
762   ierr = UpdateTS(ts,stack,e,stepnum,PETSC_FALSE);CHKERRQ(ierr);
763   PetscFunctionReturn(0);
764 }
765 
766 static PetscErrorCode SetTrajTLNR(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X)
767 {
768   Stack          *stack = &tjsch->stack;
769   PetscInt       localstepnum,laststridesize;
770   StackElement   e;
771   PetscBool      done;
772   CheckpointType cptype;
773   PetscErrorCode ierr;
774 
775   PetscFunctionBegin;
776   if (!stack->solution_only && stepnum == 0) PetscFunctionReturn(0);
777   if (stack->solution_only && stepnum == tjsch->total_steps) PetscFunctionReturn(0);
778   if (tjsch->save_stack && tjsch->recompute) PetscFunctionReturn(0);
779 
780   localstepnum = stepnum%tjsch->stride;
781   /* (stridesize-1) checkpoints are saved in each stride; an extra point is added by StackDumpAll() */
782   laststridesize = tjsch->total_steps%tjsch->stride;
783   if (!laststridesize) laststridesize = tjsch->stride;
784 
785   if (!tjsch->recompute) {
786     ierr = TopLevelStore(tj,ts,tjsch,stepnum,localstepnum,laststridesize,&done);CHKERRQ(ierr);
787     if (!tjsch->save_stack && stepnum < tjsch->total_steps-laststridesize) PetscFunctionReturn(0);
788   }
789   if (!stack->solution_only && localstepnum == 0) PetscFunctionReturn(0); /* skip last point in each stride at recompute stage or last stride */
790   if (stack->solution_only && localstepnum == tjsch->stride-1) PetscFunctionReturn(0); /* skip last step in each stride at recompute stage or last stride */
791 
792   cptype = stack->solution_only ? SOLUTIONONLY : STAGESONLY;
793   ierr = ElementCreate(ts,cptype,stack,&e);CHKERRQ(ierr);
794   ierr = ElementSet(ts,stack,&e,stepnum,time,X);CHKERRQ(ierr);
795   ierr = StackPush(stack,e);CHKERRQ(ierr);
796   PetscFunctionReturn(0);
797 }
798 
799 static PetscErrorCode GetTrajTLNR(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum)
800 {
801   Stack          *stack = &tjsch->stack;
802   PetscInt       id,localstepnum,laststridesize;
803   StackElement   e;
804   PetscErrorCode ierr;
805 
806   PetscFunctionBegin;
807   if (stepnum == tjsch->total_steps) {
808     ierr = TurnBackward(ts);CHKERRQ(ierr);
809     PetscFunctionReturn(0);
810   }
811 
812   localstepnum = stepnum%tjsch->stride;
813   laststridesize = tjsch->total_steps%tjsch->stride;
814   if (!laststridesize) laststridesize = tjsch->stride;
815   if (stack->solution_only) {
816     /* fill stack with info */
817     if (localstepnum == 0 && tjsch->total_steps-stepnum >= laststridesize) {
818       id = stepnum/tjsch->stride;
819       if (tjsch->save_stack) {
820         ierr = StackLoadAll(tj,ts,stack,id);CHKERRQ(ierr);
821         tjsch->skip_trajectory = PETSC_TRUE;
822         ierr = TurnForward(ts);CHKERRQ(ierr);
823         ierr = ReCompute(ts,tjsch,id*tjsch->stride-1,id*tjsch->stride);CHKERRQ(ierr);
824         tjsch->skip_trajectory = PETSC_FALSE;
825       } else {
826         ierr = LoadSingle(tj,ts,stack,id);CHKERRQ(ierr);
827         ierr = TurnForward(ts);CHKERRQ(ierr);
828         ierr = ReCompute(ts,tjsch,(id-1)*tjsch->stride,id*tjsch->stride);CHKERRQ(ierr);
829       }
830       PetscFunctionReturn(0);
831     }
832     /* restore a checkpoint */
833     ierr = StackPop(stack,&e);CHKERRQ(ierr);
834     ierr = UpdateTS(ts,stack,e,stepnum,PETSC_TRUE);CHKERRQ(ierr);
835     tjsch->skip_trajectory = PETSC_TRUE;
836     ierr = TurnForward(ts);CHKERRQ(ierr);
837     ierr = ReCompute(ts,tjsch,e->stepnum,stepnum);CHKERRQ(ierr);
838     tjsch->skip_trajectory = PETSC_FALSE;
839   } else {
840     CheckpointType cptype = STAGESONLY;
841     /* fill stack with info */
842     if (localstepnum == 0 && tjsch->total_steps-stepnum >= laststridesize) {
843       id = stepnum/tjsch->stride;
844       if (tjsch->save_stack) {
845         ierr = StackLoadAll(tj,ts,stack,id);CHKERRQ(ierr);
846       } else {
847         ierr = LoadSingle(tj,ts,stack,id);CHKERRQ(ierr);
848         ierr = ElementCreate(ts,cptype,stack,&e);CHKERRQ(ierr);
849         ierr = ElementSet(ts,stack,&e,(id-1)*tjsch->stride+1,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
850         ierr = StackPush(stack,e);CHKERRQ(ierr);
851         ierr = TurnForward(ts);CHKERRQ(ierr);
852         ierr = ReCompute(ts,tjsch,e->stepnum,id*tjsch->stride);CHKERRQ(ierr);
853       }
854       PetscFunctionReturn(0);
855     }
856     /* restore a checkpoint */
857     ierr = StackPop(stack,&e);CHKERRQ(ierr);
858     ierr = UpdateTS(ts,stack,e,stepnum,PETSC_TRUE);CHKERRQ(ierr);
859   }
860   PetscFunctionReturn(0);
861 }
862 
863 #if defined(PETSC_HAVE_REVOLVE)
864 static PetscErrorCode printwhattodo(PetscViewer viewer,PetscRevolveInt whattodo,RevolveCTX *rctx,PetscRevolveInt shift)
865 {
866   PetscErrorCode ierr;
867 
868   PetscFunctionBegin;
869   if (!viewer) PetscFunctionReturn(0);
870 
871   switch(whattodo) {
872     case 1:
873       ierr = PetscViewerASCIIPrintf(viewer,"Advance from %D to %D\n",rctx->oldcapo+shift,rctx->capo+shift);CHKERRQ(ierr);
874       break;
875     case 2:
876       ierr = PetscViewerASCIIPrintf(viewer,"Store in checkpoint number %D (located in RAM)\n",rctx->check);CHKERRQ(ierr);
877       break;
878     case 3:
879       ierr = PetscViewerASCIIPrintf(viewer,"First turn: Initialize adjoints and reverse first step\n");CHKERRQ(ierr);
880       break;
881     case 4:
882       ierr = PetscViewerASCIIPrintf(viewer,"Forward and reverse one step\n");CHKERRQ(ierr);
883       break;
884     case 5:
885       ierr = PetscViewerASCIIPrintf(viewer,"Restore in checkpoint number %D (located in RAM)\n",rctx->check);CHKERRQ(ierr);
886       break;
887     case 7:
888       ierr = PetscViewerASCIIPrintf(viewer,"Store in checkpoint number %D (located on disk)\n",rctx->check);CHKERRQ(ierr);
889       break;
890     case 8:
891       ierr = PetscViewerASCIIPrintf(viewer,"Restore in checkpoint number %D (located on disk)\n",rctx->check);CHKERRQ(ierr);
892       break;
893     case -1:
894       ierr = PetscViewerASCIIPrintf(viewer,"Error!");CHKERRQ(ierr);
895       break;
896   }
897   PetscFunctionReturn(0);
898 }
899 
900 static PetscErrorCode printwhattodo2(PetscViewer viewer,PetscRevolveInt whattodo,RevolveCTX *rctx,PetscRevolveInt shift)
901 {
902   PetscErrorCode ierr;
903 
904   PetscFunctionBegin;
905   if (!viewer) PetscFunctionReturn(0);
906 
907   switch(whattodo) {
908     case 1:
909       ierr = PetscViewerASCIIPrintf(viewer,"[Top Level] Advance from stride %D to stride %D\n",rctx->oldcapo+shift,rctx->capo+shift);CHKERRQ(ierr);
910       break;
911     case 2:
912       ierr = PetscViewerASCIIPrintf(viewer,"[Top Level] Store in checkpoint number %D\n",rctx->check);CHKERRQ(ierr);
913       break;
914     case 3:
915       ierr = PetscViewerASCIIPrintf(viewer,"[Top Level] First turn: Initialize adjoints and reverse first stride\n");CHKERRQ(ierr);
916       break;
917     case 4:
918       ierr = PetscViewerASCIIPrintf(viewer,"[Top Level] Forward and reverse one stride\n");CHKERRQ(ierr);
919       break;
920     case 5:
921       ierr = PetscViewerASCIIPrintf(viewer,"[Top Level] Restore in checkpoint number %D\n",rctx->check);CHKERRQ(ierr);
922       break;
923     case 7:
924       ierr = PetscViewerASCIIPrintf(viewer,"[Top Level] Store in top-level checkpoint number %D\n",rctx->check);CHKERRQ(ierr);
925       break;
926     case 8:
927       ierr = PetscViewerASCIIPrintf(viewer,"[Top Level] Restore in top-level checkpoint number %D\n",rctx->check);CHKERRQ(ierr);
928       break;
929     case -1:
930       ierr = PetscViewerASCIIPrintf(viewer,"[Top Level] Error!");CHKERRQ(ierr);
931       break;
932   }
933   PetscFunctionReturn(0);
934 }
935 
936 static PetscErrorCode InitRevolve(PetscInt fine,PetscInt snaps,RevolveCTX *rctx)
937 {
938   PetscRevolveInt rsnaps,rfine;
939   PetscErrorCode  ierr;
940 
941   PetscFunctionBegin;
942   ierr = PetscRevolveIntCast(snaps,&rsnaps);CHKERRQ(ierr);
943   ierr = PetscRevolveIntCast(fine,&rfine);CHKERRQ(ierr);
944   revolve_reset();
945   revolve_create_offline(rfine,rsnaps);
946   rctx->snaps_in       = rsnaps;
947   rctx->fine           = rfine;
948   rctx->check          = 0;
949   rctx->capo           = 0;
950   rctx->reverseonestep = PETSC_FALSE;
951   /* check stepsleft? */
952   PetscFunctionReturn(0);
953 }
954 
955 static PetscErrorCode FastForwardRevolve(RevolveCTX *rctx)
956 {
957   PetscRevolveInt whattodo;
958 
959   PetscFunctionBegin;
960   whattodo = 0;
961   while (whattodo!=3) { /* we have to fast forward revolve to the beginning of the backward sweep due to unfriendly revolve interface */
962     whattodo = revolve_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where);
963   }
964   PetscFunctionReturn(0);
965 }
966 
967 static PetscErrorCode ApplyRevolve(PetscViewer viewer,SchedulerType stype,RevolveCTX *rctx,PetscRevolveInt total_steps,PetscRevolveInt stepnum,PetscRevolveInt localstepnum,PetscBool toplevel,PetscInt *store)
968 {
969   PetscErrorCode  ierr;
970   PetscRevolveInt shift,whattodo;
971 
972   PetscFunctionBegin;
973   *store = 0;
974   if (rctx->stepsleft > 0) { /* advance the solution without checkpointing anything as Revolve requires */
975     rctx->stepsleft--;
976     PetscFunctionReturn(0);
977   }
978   /* let Revolve determine what to do next */
979   shift         = stepnum-localstepnum;
980   rctx->oldcapo = rctx->capo;
981   rctx->capo    = localstepnum;
982 
983   if (!toplevel) whattodo = revolve_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where);
984   else whattodo = revolve2_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where);
985   if (stype == REVOLVE_ONLINE && whattodo == 8) whattodo = 5;
986   if (stype == REVOLVE_ONLINE && whattodo == 7) whattodo = 2;
987   if (!toplevel) {ierr = printwhattodo(viewer,whattodo,rctx,shift);CHKERRQ(ierr);}
988   else {ierr = printwhattodo2(viewer,whattodo,rctx,shift);CHKERRQ(ierr);}
989   if (whattodo == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in the Revolve library");
990   if (whattodo == 1) { /* advance some time steps */
991     if (stype == REVOLVE_ONLINE && rctx->capo >= total_steps-1) {
992       revolve_turn(total_steps,&rctx->capo,&rctx->fine);
993       if (!toplevel) {ierr = printwhattodo(viewer,whattodo,rctx,shift);CHKERRQ(ierr);}
994       else {ierr = printwhattodo2(viewer,whattodo,rctx,shift);CHKERRQ(ierr);}
995     }
996     rctx->stepsleft = rctx->capo-rctx->oldcapo-1;
997   }
998   if (whattodo == 3 || whattodo == 4) { /* ready for a reverse step */
999     rctx->reverseonestep = PETSC_TRUE;
1000   }
1001   if (whattodo == 5) { /* restore a checkpoint and ask Revolve what to do next */
1002     rctx->oldcapo = rctx->capo;
1003     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*/
1004     else whattodo = revolve2_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where);
1005     if (!toplevel) {ierr = printwhattodo(viewer,whattodo,rctx,shift);CHKERRQ(ierr);}
1006     else {ierr = printwhattodo2(viewer,whattodo,rctx,shift);CHKERRQ(ierr);}
1007     if (whattodo == 3 || whattodo == 4) rctx->reverseonestep = PETSC_TRUE;
1008     if (whattodo == 1) rctx->stepsleft = rctx->capo-rctx->oldcapo;
1009   }
1010   if (whattodo == 7) { /* save the checkpoint to disk */
1011     *store = 2;
1012     rctx->oldcapo = rctx->capo;
1013     whattodo = revolve_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where); /* must return 1 */
1014     ierr = printwhattodo(viewer,whattodo,rctx,shift);CHKERRQ(ierr);
1015     rctx->stepsleft = rctx->capo-rctx->oldcapo-1;
1016   }
1017   if (whattodo == 2) { /* store a checkpoint to RAM and ask Revolve how many time steps to advance next */
1018     *store = 1;
1019     rctx->oldcapo = rctx->capo;
1020     if (!toplevel) whattodo = revolve_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where); /* must return 1 */
1021     else whattodo = revolve2_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where);
1022     if (!toplevel) {ierr = printwhattodo(viewer,whattodo,rctx,shift);CHKERRQ(ierr);}
1023     else {ierr = printwhattodo2(viewer,whattodo,rctx,shift);CHKERRQ(ierr);}
1024     if (stype == REVOLVE_ONLINE && rctx->capo >= total_steps-1) {
1025       revolve_turn(total_steps,&rctx->capo,&rctx->fine);
1026       ierr = printwhattodo(viewer,whattodo,rctx,shift);CHKERRQ(ierr);
1027     }
1028     rctx->stepsleft = rctx->capo-rctx->oldcapo-1;
1029   }
1030   PetscFunctionReturn(0);
1031 }
1032 
1033 static PetscErrorCode SetTrajROF(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X)
1034 {
1035   Stack           *stack = &tjsch->stack;
1036   PetscInt        store;
1037   StackElement    e;
1038   PetscRevolveInt rtotal_steps,rstepnum;
1039   CheckpointType  cptype;
1040   PetscErrorCode  ierr;
1041 
1042   PetscFunctionBegin;
1043   if (!stack->solution_only && stepnum == 0) PetscFunctionReturn(0);
1044   if (stack->solution_only && stepnum == tjsch->total_steps) PetscFunctionReturn(0);
1045   ierr = PetscRevolveIntCast(tjsch->total_steps,&rtotal_steps);CHKERRQ(ierr);
1046   ierr = PetscRevolveIntCast(stepnum,&rstepnum);CHKERRQ(ierr);
1047   ierr = ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,rtotal_steps,rstepnum,rstepnum,PETSC_FALSE,&store);CHKERRQ(ierr);
1048   if (store == 1) {
1049     if (stepnum < stack->top) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_MEMC,"Illegal modification of a non-top stack element");
1050     cptype = stack->solution_only ? SOLUTIONONLY : SOLUTION_STAGES;
1051     ierr = ElementCreate(ts,cptype,stack,&e);CHKERRQ(ierr);
1052     ierr = ElementSet(ts,stack,&e,stepnum,time,X);CHKERRQ(ierr);
1053     ierr = StackPush(stack,e);CHKERRQ(ierr);
1054   }
1055   PetscFunctionReturn(0);
1056 }
1057 
1058 static PetscErrorCode GetTrajROF(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum)
1059 {
1060   Stack           *stack = &tjsch->stack;
1061   PetscInt        store;
1062   PetscRevolveInt whattodo,shift,rtotal_steps,rstepnum;
1063   StackElement    e;
1064   PetscErrorCode  ierr;
1065 
1066   PetscFunctionBegin;
1067   if (stepnum == 0 || stepnum == tjsch->total_steps) {
1068     ierr = TurnBackward(ts);CHKERRQ(ierr);
1069     tjsch->rctx->reverseonestep = PETSC_FALSE;
1070     PetscFunctionReturn(0);
1071   }
1072   /* restore a checkpoint */
1073   ierr = StackTop(stack,&e);CHKERRQ(ierr);
1074   ierr = UpdateTS(ts,stack,e,stepnum,PETSC_TRUE);CHKERRQ(ierr);
1075   ierr = PetscRevolveIntCast(tjsch->total_steps,&rtotal_steps);CHKERRQ(ierr);
1076   ierr = PetscRevolveIntCast(stepnum,&rstepnum);CHKERRQ(ierr);
1077   if (stack->solution_only) { /* start with restoring a checkpoint */
1078     tjsch->rctx->capo = rstepnum;
1079     tjsch->rctx->oldcapo = tjsch->rctx->capo;
1080     shift = 0;
1081     whattodo = revolve_action(&tjsch->rctx->check,&tjsch->rctx->capo,&tjsch->rctx->fine,tjsch->rctx->snaps_in,&tjsch->rctx->info,&tjsch->rctx->where);
1082     ierr = printwhattodo(tj->monitor,whattodo,tjsch->rctx,shift);CHKERRQ(ierr);
1083   } else { /* 2 revolve actions: restore a checkpoint and then advance */
1084     ierr = ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,rtotal_steps,rstepnum,rstepnum,PETSC_FALSE,&store);CHKERRQ(ierr);
1085     if (tj->monitor) {
1086       ierr = PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr);
1087       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);
1088       ierr = PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr);
1089     }
1090     if (!tjsch->rctx->reverseonestep && tjsch->rctx->stepsleft > 0) tjsch->rctx->stepsleft--;
1091   }
1092   if (stack->solution_only || (!stack->solution_only && e->stepnum < stepnum)) {
1093     ierr = TurnForward(ts);CHKERRQ(ierr);
1094     ierr = ReCompute(ts,tjsch,e->stepnum,stepnum);CHKERRQ(ierr);
1095   }
1096   if ((stack->solution_only && e->stepnum+1 == stepnum) || (!stack->solution_only && e->stepnum == stepnum)) {
1097     ierr = StackPop(stack,&e);CHKERRQ(ierr);
1098   }
1099   tjsch->rctx->reverseonestep = PETSC_FALSE;
1100   PetscFunctionReturn(0);
1101 }
1102 
1103 static PetscErrorCode SetTrajRON(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X)
1104 {
1105   Stack           *stack = &tjsch->stack;
1106   Vec             *Y;
1107   PetscInt        i,store;
1108   PetscReal       timeprev;
1109   StackElement    e;
1110   RevolveCTX      *rctx = tjsch->rctx;
1111   PetscRevolveInt rtotal_steps,rstepnum;
1112   CheckpointType  cptype;
1113   PetscErrorCode  ierr;
1114 
1115   PetscFunctionBegin;
1116   if (!stack->solution_only && stepnum == 0) PetscFunctionReturn(0);
1117   if (stack->solution_only && stepnum == tjsch->total_steps) PetscFunctionReturn(0);
1118   ierr = PetscRevolveIntCast(tjsch->total_steps,&rtotal_steps);CHKERRQ(ierr);
1119   ierr = PetscRevolveIntCast(stepnum,&rstepnum);CHKERRQ(ierr);
1120   ierr = ApplyRevolve(tj->monitor,tjsch->stype,rctx,rtotal_steps,rstepnum,rstepnum,PETSC_FALSE,&store);CHKERRQ(ierr);
1121   if (store == 1) {
1122     if (rctx->check != stack->top+1) { /* overwrite some non-top checkpoint in the stack */
1123       ierr = StackFind(stack,&e,rctx->check);CHKERRQ(ierr);
1124       if (HaveSolution(e->cptype)) {
1125         ierr = VecCopy(X,e->X);CHKERRQ(ierr);
1126       }
1127       if (HaveStages(e->cptype)) {
1128         ierr = TSGetStages(ts,&stack->numY,&Y);CHKERRQ(ierr);
1129         for (i=0;i<stack->numY;i++) {
1130           ierr = VecCopy(Y[i],e->Y[i]);CHKERRQ(ierr);
1131         }
1132       }
1133       e->stepnum  = stepnum;
1134       e->time     = time;
1135       ierr        = TSGetPrevTime(ts,&timeprev);CHKERRQ(ierr);
1136       e->timeprev = timeprev;
1137     } else {
1138       if (stepnum < stack->top) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_MEMC,"Illegal modification of a non-top stack element");
1139       cptype = stack->solution_only ? SOLUTIONONLY : SOLUTION_STAGES;
1140       ierr = ElementCreate(ts,cptype,stack,&e);CHKERRQ(ierr);
1141       ierr = ElementSet(ts,stack,&e,stepnum,time,X);CHKERRQ(ierr);
1142       ierr = StackPush(stack,e);CHKERRQ(ierr);
1143     }
1144   }
1145   PetscFunctionReturn(0);
1146 }
1147 
1148 static PetscErrorCode GetTrajRON(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum)
1149 {
1150   Stack           *stack = &tjsch->stack;
1151   PetscRevolveInt whattodo,shift,rstepnum;
1152   StackElement    e;
1153   PetscErrorCode  ierr;
1154 
1155   PetscFunctionBegin;
1156   if (stepnum == 0 || stepnum == tjsch->total_steps) {
1157     ierr = TurnBackward(ts);CHKERRQ(ierr);
1158     tjsch->rctx->reverseonestep = PETSC_FALSE;
1159     PetscFunctionReturn(0);
1160   }
1161   ierr = PetscRevolveIntCast(stepnum,&rstepnum);CHKERRQ(ierr);
1162   tjsch->rctx->capo = rstepnum;
1163   tjsch->rctx->oldcapo = tjsch->rctx->capo;
1164   shift = 0;
1165   whattodo = revolve_action(&tjsch->rctx->check,&tjsch->rctx->capo,&tjsch->rctx->fine,tjsch->rctx->snaps_in,&tjsch->rctx->info,&tjsch->rctx->where); /* whattodo=restore */
1166   if (whattodo == 8) whattodo = 5;
1167   ierr = printwhattodo(tj->monitor,whattodo,tjsch->rctx,shift);CHKERRQ(ierr);
1168   /* restore a checkpoint */
1169   ierr = StackFind(stack,&e,tjsch->rctx->check);CHKERRQ(ierr);
1170   ierr = UpdateTS(ts,stack,e,stepnum,PETSC_TRUE);CHKERRQ(ierr);
1171   if (!stack->solution_only) { /* whattodo must be 5 */
1172     /* ask Revolve what to do next */
1173     tjsch->rctx->oldcapo = tjsch->rctx->capo;
1174     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*/
1175     ierr = printwhattodo(tj->monitor,whattodo,tjsch->rctx,shift);CHKERRQ(ierr);
1176     if (whattodo == 3 || whattodo == 4) tjsch->rctx->reverseonestep = PETSC_TRUE;
1177     if (whattodo == 1) tjsch->rctx->stepsleft = tjsch->rctx->capo-tjsch->rctx->oldcapo;
1178     if (tj->monitor) {
1179       ierr = PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr);
1180       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);
1181       ierr = PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr);
1182     }
1183     if (!tjsch->rctx->reverseonestep && tjsch->rctx->stepsleft > 0) tjsch->rctx->stepsleft--;
1184   }
1185   if (stack->solution_only || (!stack->solution_only && e->stepnum < stepnum)) {
1186     ierr = TurnForward(ts);CHKERRQ(ierr);
1187     ierr = ReCompute(ts,tjsch,e->stepnum,stepnum);CHKERRQ(ierr);
1188   }
1189   tjsch->rctx->reverseonestep = PETSC_FALSE;
1190   PetscFunctionReturn(0);
1191 }
1192 
1193 static PetscErrorCode SetTrajTLR(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X)
1194 {
1195   Stack           *stack = &tjsch->stack;
1196   PetscInt        store,localstepnum,laststridesize;
1197   StackElement    e;
1198   PetscBool       done = PETSC_FALSE;
1199   PetscRevolveInt rtotal_steps,rstepnum,rlocalstepnum;
1200   CheckpointType  cptype;
1201   PetscErrorCode  ierr;
1202 
1203   PetscFunctionBegin;
1204   if (!stack->solution_only && stepnum == 0) PetscFunctionReturn(0);
1205   if (stack->solution_only && stepnum == tjsch->total_steps) PetscFunctionReturn(0);
1206 
1207   localstepnum = stepnum%tjsch->stride;
1208   laststridesize = tjsch->total_steps%tjsch->stride;
1209   if (!laststridesize) laststridesize = tjsch->stride;
1210 
1211   if (!tjsch->recompute) {
1212     ierr = TopLevelStore(tj,ts,tjsch,stepnum,localstepnum,laststridesize,&done);CHKERRQ(ierr);
1213     /* revolve is needed for the last stride; different starting points for last stride between solutin_only and !solutin_only */
1214     if (!stack->solution_only && !tjsch->save_stack && stepnum <= tjsch->total_steps-laststridesize) PetscFunctionReturn(0);
1215     if (stack->solution_only && !tjsch->save_stack && stepnum < tjsch->total_steps-laststridesize) PetscFunctionReturn(0);
1216   }
1217   if (tjsch->save_stack && done) {
1218     ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr);
1219     PetscFunctionReturn(0);
1220   }
1221   if (laststridesize < tjsch->stride) {
1222     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 */
1223       ierr = InitRevolve(laststridesize,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr);
1224     }
1225     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 */
1226       ierr = InitRevolve(laststridesize,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr);
1227     }
1228   }
1229   ierr = PetscRevolveIntCast(tjsch->total_steps,&rtotal_steps);CHKERRQ(ierr);
1230   ierr = PetscRevolveIntCast(stepnum,&rstepnum);CHKERRQ(ierr);
1231   ierr = PetscRevolveIntCast(localstepnum,&rlocalstepnum);CHKERRQ(ierr);
1232   ierr = ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,rtotal_steps,rstepnum,rlocalstepnum,PETSC_FALSE,&store);CHKERRQ(ierr);
1233   if (store == 1) {
1234     if (localstepnum < stack->top) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_MEMC,"Illegal modification of a non-top stack element");
1235     cptype = stack->solution_only ? SOLUTIONONLY : SOLUTION_STAGES;
1236     ierr = ElementCreate(ts,cptype,stack,&e);CHKERRQ(ierr);
1237     ierr = ElementSet(ts,stack,&e,stepnum,time,X);CHKERRQ(ierr);
1238     ierr = StackPush(stack,e);CHKERRQ(ierr);
1239   }
1240   PetscFunctionReturn(0);
1241 }
1242 
1243 static PetscErrorCode GetTrajTLR(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum)
1244 {
1245   Stack           *stack = &tjsch->stack;
1246   PetscRevolveInt whattodo,shift,rstepnum,rlocalstepnum,rtotal_steps;
1247   PetscInt        localstepnum,stridenum,laststridesize,store;
1248   StackElement    e;
1249   CheckpointType  cptype;
1250   PetscErrorCode  ierr;
1251 
1252   PetscFunctionBegin;
1253   localstepnum = stepnum%tjsch->stride;
1254   stridenum    = stepnum/tjsch->stride;
1255   if (stepnum == tjsch->total_steps) {
1256     ierr = TurnBackward(ts);CHKERRQ(ierr);
1257     tjsch->rctx->reverseonestep = PETSC_FALSE;
1258     PetscFunctionReturn(0);
1259   }
1260   laststridesize = tjsch->total_steps%tjsch->stride;
1261   if (!laststridesize) laststridesize = tjsch->stride;
1262   ierr = PetscRevolveIntCast(tjsch->total_steps,&rtotal_steps);CHKERRQ(ierr);
1263   ierr = PetscRevolveIntCast(stepnum,&rstepnum);CHKERRQ(ierr);
1264   ierr = PetscRevolveIntCast(localstepnum,&rlocalstepnum);CHKERRQ(ierr);
1265   if (stack->solution_only) {
1266     /* fill stack */
1267     if (localstepnum == 0 && stepnum <= tjsch->total_steps-laststridesize) {
1268       if (tjsch->save_stack) {
1269         ierr = StackLoadAll(tj,ts,stack,stridenum);CHKERRQ(ierr);
1270         ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr);
1271         ierr = FastForwardRevolve(tjsch->rctx);CHKERRQ(ierr);
1272         tjsch->skip_trajectory = PETSC_TRUE;
1273         ierr = TurnForward(ts);CHKERRQ(ierr);
1274         ierr = ReCompute(ts,tjsch,stridenum*tjsch->stride-1,stridenum*tjsch->stride);CHKERRQ(ierr);
1275         tjsch->skip_trajectory = PETSC_FALSE;
1276       } else {
1277         ierr = LoadSingle(tj,ts,stack,stridenum);CHKERRQ(ierr);
1278         ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr);
1279         ierr = TurnForward(ts);CHKERRQ(ierr);
1280         ierr = ReCompute(ts,tjsch,(stridenum-1)*tjsch->stride,stridenum*tjsch->stride);CHKERRQ(ierr);
1281       }
1282       PetscFunctionReturn(0);
1283     }
1284     /* restore a checkpoint */
1285     ierr = StackTop(stack,&e);CHKERRQ(ierr);
1286     ierr = UpdateTS(ts,stack,e,stepnum,PETSC_TRUE);CHKERRQ(ierr);
1287     /* start with restoring a checkpoint */
1288     tjsch->rctx->capo = rstepnum;
1289     tjsch->rctx->oldcapo = tjsch->rctx->capo;
1290     shift = rstepnum-rlocalstepnum;
1291     whattodo = revolve_action(&tjsch->rctx->check,&tjsch->rctx->capo,&tjsch->rctx->fine,tjsch->rctx->snaps_in,&tjsch->rctx->info,&tjsch->rctx->where);
1292     ierr = printwhattodo(tj->monitor,whattodo,tjsch->rctx,shift);CHKERRQ(ierr);
1293     ierr = TurnForward(ts);CHKERRQ(ierr);
1294     ierr = ReCompute(ts,tjsch,e->stepnum,stepnum);CHKERRQ(ierr);
1295     if (e->stepnum+1 == stepnum) {
1296       ierr = StackPop(stack,&e);CHKERRQ(ierr);
1297     }
1298   } else {
1299     /* fill stack with info */
1300     if (localstepnum == 0 && tjsch->total_steps-stepnum >= laststridesize) {
1301       if (tjsch->save_stack) {
1302         ierr = StackLoadAll(tj,ts,stack,stridenum);CHKERRQ(ierr);
1303         ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr);
1304         ierr = FastForwardRevolve(tjsch->rctx);CHKERRQ(ierr);
1305       } else {
1306         PetscRevolveInt rnum;
1307         ierr = LoadSingle(tj,ts,stack,stridenum);CHKERRQ(ierr);
1308         ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr);
1309         ierr = PetscRevolveIntCast((stridenum-1)*tjsch->stride+1,&rnum);CHKERRQ(ierr);
1310         ierr = ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,rtotal_steps,rnum,1,PETSC_FALSE,&store);CHKERRQ(ierr);
1311         if (tj->monitor) {
1312           ierr = PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr);
1313           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);
1314           ierr = PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr);
1315         }
1316         cptype = SOLUTION_STAGES;
1317         ierr = ElementCreate(ts,cptype,stack,&e);CHKERRQ(ierr);
1318         ierr = ElementSet(ts,stack,&e,(stridenum-1)*tjsch->stride+1,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
1319         ierr = StackPush(stack,e);CHKERRQ(ierr);
1320         ierr = TurnForward(ts);CHKERRQ(ierr);
1321         ierr = ReCompute(ts,tjsch,e->stepnum,stridenum*tjsch->stride);CHKERRQ(ierr);
1322       }
1323       PetscFunctionReturn(0);
1324     }
1325     /* restore a checkpoint */
1326     ierr = StackTop(stack,&e);CHKERRQ(ierr);
1327     ierr = UpdateTS(ts,stack,e,stepnum,PETSC_TRUE);CHKERRQ(ierr);
1328     /* 2 revolve actions: restore a checkpoint and then advance */
1329     ierr = ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,rtotal_steps,rstepnum,rlocalstepnum,PETSC_FALSE,&store);CHKERRQ(ierr);
1330     if (tj->monitor) {
1331       ierr = PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr);
1332       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);
1333       ierr = PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr);
1334     }
1335     if (!tjsch->rctx->reverseonestep && tjsch->rctx->stepsleft > 0) tjsch->rctx->stepsleft--;
1336     if (e->stepnum < stepnum) {
1337       ierr = TurnForward(ts);CHKERRQ(ierr);
1338       ierr = ReCompute(ts,tjsch,e->stepnum,stepnum);CHKERRQ(ierr);
1339     }
1340     if (e->stepnum == stepnum) {
1341       ierr = StackPop(stack,&e);CHKERRQ(ierr);
1342     }
1343   }
1344   tjsch->rctx->reverseonestep = PETSC_FALSE;
1345   PetscFunctionReturn(0);
1346 }
1347 
1348 static PetscErrorCode SetTrajTLTR(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X)
1349 {
1350   Stack           *stack = &tjsch->stack;
1351   PetscInt        store,localstepnum,stridenum,laststridesize;
1352   StackElement    e;
1353   PetscBool       done = PETSC_FALSE;
1354   PetscRevolveInt rlocalstepnum,rstepnum,rtotal_steps;
1355   PetscErrorCode  ierr;
1356 
1357   PetscFunctionBegin;
1358   if (!stack->solution_only && stepnum == 0) PetscFunctionReturn(0);
1359   if (stack->solution_only && stepnum == tjsch->total_steps) PetscFunctionReturn(0);
1360 
1361   localstepnum = stepnum%tjsch->stride; /* index at the bottom level (inside a stride) */
1362   stridenum    = stepnum/tjsch->stride; /* index at the top level */
1363   laststridesize = tjsch->total_steps%tjsch->stride;
1364   if (!laststridesize) laststridesize = tjsch->stride;
1365   if (stack->solution_only && localstepnum == 0 && !tjsch->rctx2->reverseonestep) {
1366     ierr = PetscRevolveIntCast((tjsch->total_steps+tjsch->stride-1)/tjsch->stride,&rtotal_steps);CHKERRQ(ierr);
1367     ierr = PetscRevolveIntCast(stridenum,&rstepnum);CHKERRQ(ierr);
1368     ierr = ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx2,rtotal_steps,rstepnum,rstepnum,PETSC_TRUE,&tjsch->store_stride);CHKERRQ(ierr);
1369     if (laststridesize < tjsch->stride && stepnum == tjsch->total_steps-laststridesize) {
1370       ierr = InitRevolve(laststridesize,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr);
1371     }
1372   }
1373   if (!stack->solution_only && localstepnum == 1 && !tjsch->rctx2->reverseonestep) {
1374     ierr = PetscRevolveIntCast((tjsch->total_steps+tjsch->stride-1)/tjsch->stride,&rtotal_steps);CHKERRQ(ierr);
1375     ierr = PetscRevolveIntCast(stridenum,&rstepnum);CHKERRQ(ierr);
1376     ierr = ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx2,rtotal_steps,rstepnum,rstepnum,PETSC_TRUE,&tjsch->store_stride);CHKERRQ(ierr);
1377     if (laststridesize < tjsch->stride && stepnum == tjsch->total_steps-laststridesize+1) {
1378       ierr = InitRevolve(laststridesize,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr);
1379     }
1380   }
1381   if (tjsch->store_stride) {
1382     ierr = TopLevelStore(tj,ts,tjsch,stepnum,localstepnum,laststridesize,&done);CHKERRQ(ierr);
1383     if (done) {
1384       ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr);
1385       PetscFunctionReturn(0);
1386     }
1387   }
1388   if (stepnum < tjsch->total_steps-laststridesize) {
1389     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 */
1390     if (!tjsch->save_stack && !tjsch->rctx2->reverseonestep) PetscFunctionReturn(0); /* store operation does not require revolve be called at bottom level */
1391   }
1392   /* 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 */
1393   if (!stack->solution_only && localstepnum == 0 && stepnum != tjsch->total_steps && !tjsch->recompute) PetscFunctionReturn(0);
1394   ierr = PetscRevolveIntCast(tjsch->total_steps,&rtotal_steps);CHKERRQ(ierr);
1395   ierr = PetscRevolveIntCast(stepnum,&rstepnum);CHKERRQ(ierr);
1396   ierr = PetscRevolveIntCast(localstepnum,&rlocalstepnum);CHKERRQ(ierr);
1397   ierr = ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,rtotal_steps,rstepnum,rlocalstepnum,PETSC_FALSE,&store);CHKERRQ(ierr);
1398   if (store == 1) {
1399     CheckpointType cptype;
1400     if (localstepnum < stack->top) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_MEMC,"Illegal modification of a non-top stack element");
1401     cptype = stack->solution_only ? SOLUTIONONLY : SOLUTION_STAGES;
1402     ierr = ElementCreate(ts,cptype,stack,&e);CHKERRQ(ierr);
1403     ierr = ElementSet(ts,stack,&e,stepnum,time,X);CHKERRQ(ierr);
1404     ierr = StackPush(stack,e);CHKERRQ(ierr);
1405   }
1406   PetscFunctionReturn(0);
1407 }
1408 
1409 static PetscErrorCode GetTrajTLTR(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum)
1410 {
1411   Stack           *stack = &tjsch->stack;
1412   DiskStack       *diskstack = &tjsch->diskstack;
1413   PetscInt        localstepnum,stridenum,restoredstridenum,laststridesize,store;
1414   StackElement    e;
1415   PetscRevolveInt whattodo,shift;
1416   PetscRevolveInt rtotal_steps,rstepnum,rlocalstepnum;
1417   PetscErrorCode  ierr;
1418 
1419   PetscFunctionBegin;
1420   localstepnum = stepnum%tjsch->stride;
1421   stridenum    = stepnum/tjsch->stride;
1422   if (stepnum == tjsch->total_steps) {
1423     ierr = TurnBackward(ts);CHKERRQ(ierr);
1424     tjsch->rctx->reverseonestep = PETSC_FALSE;
1425     PetscFunctionReturn(0);
1426   }
1427   laststridesize = tjsch->total_steps%tjsch->stride;
1428   if (!laststridesize) laststridesize = tjsch->stride;
1429   /*
1430    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:
1431      Case 1 (save_stack)
1432        Restore a disk checkpoint; update TS with the last element in the restored data; recompute to the current point.
1433      Case 2 (!save_stack)
1434        Restore a disk checkpoint; update TS with the restored point; recompute to the current point.
1435   */
1436   if (localstepnum == 0 && stepnum <= tjsch->total_steps-laststridesize) {
1437     /* restore the top element in the stack for disk checkpoints */
1438     restoredstridenum = diskstack->container[diskstack->top];
1439     tjsch->rctx2->reverseonestep = PETSC_FALSE;
1440     /* top-level revolve must be applied before current step, just like the solution_only mode for single-level revolve */
1441     if (!tjsch->save_stack && stack->solution_only) { /* start with restoring a checkpoint */
1442       ierr = PetscRevolveIntCast(stridenum,&rstepnum);CHKERRQ(ierr);
1443       tjsch->rctx2->capo = rstepnum;
1444       tjsch->rctx2->oldcapo = tjsch->rctx2->capo;
1445       shift = 0;
1446       whattodo = revolve2_action(&tjsch->rctx2->check,&tjsch->rctx2->capo,&tjsch->rctx2->fine,tjsch->rctx2->snaps_in,&tjsch->rctx2->info,&tjsch->rctx2->where);
1447       ierr = printwhattodo2(tj->monitor,whattodo,tjsch->rctx2,shift);CHKERRQ(ierr);
1448     } else { /* 2 revolve actions: restore a checkpoint and then advance */
1449       ierr = PetscRevolveIntCast((tjsch->total_steps+tjsch->stride-1)/tjsch->stride,&rtotal_steps);CHKERRQ(ierr);
1450       ierr = PetscRevolveIntCast(stridenum,&rstepnum);CHKERRQ(ierr);
1451       ierr = ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx2,rtotal_steps,rstepnum,rstepnum,PETSC_TRUE,&tjsch->store_stride);CHKERRQ(ierr);
1452       if (tj->monitor) {
1453         ierr = PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr);
1454         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);
1455         ierr = PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr);
1456       }
1457       if (!tjsch->rctx2->reverseonestep && tjsch->rctx2->stepsleft > 0) tjsch->rctx2->stepsleft--;
1458     }
1459     /* fill stack */
1460     if (stack->solution_only) {
1461       if (tjsch->save_stack) {
1462         if (restoredstridenum < stridenum) {
1463           ierr = StackLoadLast(tj,ts,stack,restoredstridenum);CHKERRQ(ierr);
1464         } else {
1465           ierr = StackLoadAll(tj,ts,stack,restoredstridenum);CHKERRQ(ierr);
1466         }
1467         /* recompute one step ahead */
1468         tjsch->skip_trajectory = PETSC_TRUE;
1469         ierr = TurnForward(ts);CHKERRQ(ierr);
1470         ierr = ReCompute(ts,tjsch,stridenum*tjsch->stride-1,stridenum*tjsch->stride);CHKERRQ(ierr);
1471         tjsch->skip_trajectory = PETSC_FALSE;
1472         if (restoredstridenum < stridenum) {
1473           ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr);
1474           ierr = TurnForward(ts);CHKERRQ(ierr);
1475           ierr = ReCompute(ts,tjsch,restoredstridenum*tjsch->stride,stepnum);CHKERRQ(ierr);
1476         } else { /* stack ready, fast forward revolve status */
1477           ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr);
1478           ierr = FastForwardRevolve(tjsch->rctx);CHKERRQ(ierr);
1479         }
1480       } else {
1481         ierr = LoadSingle(tj,ts,stack,restoredstridenum);CHKERRQ(ierr);
1482         ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr);
1483         ierr = TurnForward(ts);CHKERRQ(ierr);
1484         ierr = ReCompute(ts,tjsch,(restoredstridenum-1)*tjsch->stride,stepnum);CHKERRQ(ierr);
1485       }
1486     } else {
1487       if (tjsch->save_stack) {
1488         if (restoredstridenum < stridenum) {
1489           ierr = StackLoadLast(tj,ts,stack,restoredstridenum);CHKERRQ(ierr);
1490           /* reset revolve */
1491           ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr);
1492           ierr = TurnForward(ts);CHKERRQ(ierr);
1493           ierr = ReCompute(ts,tjsch,restoredstridenum*tjsch->stride,stepnum);CHKERRQ(ierr);
1494         } else { /* stack ready, fast forward revolve status */
1495           ierr = StackLoadAll(tj,ts,stack,restoredstridenum);CHKERRQ(ierr);
1496           ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr);
1497           ierr = FastForwardRevolve(tjsch->rctx);CHKERRQ(ierr);
1498         }
1499       } else {
1500         ierr = LoadSingle(tj,ts,stack,restoredstridenum);CHKERRQ(ierr);
1501         ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr);
1502         /* push first element to stack */
1503         if (tjsch->store_stride || tjsch->rctx2->reverseonestep) {
1504           CheckpointType cptype = SOLUTION_STAGES;
1505           shift = (restoredstridenum-1)*tjsch->stride-localstepnum;
1506           ierr = PetscRevolveIntCast(tjsch->total_steps,&rtotal_steps);CHKERRQ(ierr);
1507           ierr = PetscRevolveIntCast((restoredstridenum-1)*tjsch->stride+1,&rstepnum);CHKERRQ(ierr);
1508           ierr = ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,rtotal_steps,rstepnum,1,PETSC_FALSE,&store);CHKERRQ(ierr);
1509           if (tj->monitor) {
1510             ierr = PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr);
1511             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);
1512             ierr = PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr);
1513           }
1514           ierr = ElementCreate(ts,cptype,stack,&e);CHKERRQ(ierr);
1515           ierr = ElementSet(ts,stack,&e,(restoredstridenum-1)*tjsch->stride+1,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
1516           ierr = StackPush(stack,e);CHKERRQ(ierr);
1517         }
1518         ierr = TurnForward(ts);CHKERRQ(ierr);
1519         ierr = ReCompute(ts,tjsch,(restoredstridenum-1)*tjsch->stride+1,stepnum);CHKERRQ(ierr);
1520       }
1521     }
1522     if (restoredstridenum == stridenum) diskstack->top--;
1523     tjsch->rctx->reverseonestep = PETSC_FALSE;
1524     PetscFunctionReturn(0);
1525   }
1526 
1527   if (stack->solution_only) {
1528     /* restore a checkpoint */
1529     ierr = StackTop(stack,&e);CHKERRQ(ierr);
1530     ierr = UpdateTS(ts,stack,e,stepnum,PETSC_TRUE);CHKERRQ(ierr);
1531     /* start with restoring a checkpoint */
1532     ierr = PetscRevolveIntCast(stepnum,&rstepnum);CHKERRQ(ierr);
1533     ierr = PetscRevolveIntCast(localstepnum,&rlocalstepnum);CHKERRQ(ierr);
1534     tjsch->rctx->capo = rstepnum;
1535     tjsch->rctx->oldcapo = tjsch->rctx->capo;
1536     shift = rstepnum-rlocalstepnum;
1537     whattodo = revolve_action(&tjsch->rctx->check,&tjsch->rctx->capo,&tjsch->rctx->fine,tjsch->rctx->snaps_in,&tjsch->rctx->info,&tjsch->rctx->where);
1538     ierr = printwhattodo(tj->monitor,whattodo,tjsch->rctx,shift);CHKERRQ(ierr);
1539     ierr = TurnForward(ts);CHKERRQ(ierr);
1540     ierr = ReCompute(ts,tjsch,e->stepnum,stepnum);CHKERRQ(ierr);
1541     if (e->stepnum+1 == stepnum) {
1542       ierr = StackPop(stack,&e);CHKERRQ(ierr);
1543     }
1544   } else {
1545     PetscRevolveInt rlocalstepnum;
1546     /* restore a checkpoint */
1547     ierr = StackTop(stack,&e);CHKERRQ(ierr);
1548     ierr = UpdateTS(ts,stack,e,stepnum,PETSC_TRUE);CHKERRQ(ierr);
1549     /* 2 revolve actions: restore a checkpoint and then advance */
1550     ierr = PetscRevolveIntCast(tjsch->total_steps,&rtotal_steps);CHKERRQ(ierr);
1551     ierr = PetscRevolveIntCast(stridenum,&rstepnum);CHKERRQ(ierr);
1552     ierr = PetscRevolveIntCast(localstepnum,&rlocalstepnum);CHKERRQ(ierr);
1553     ierr = ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,rtotal_steps,rstepnum,rlocalstepnum,PETSC_FALSE,&store);CHKERRQ(ierr);
1554     if (tj->monitor) {
1555       ierr = PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr);
1556       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);
1557       ierr = PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr);
1558     }
1559     if (!tjsch->rctx->reverseonestep && tjsch->rctx->stepsleft > 0) tjsch->rctx->stepsleft--;
1560     if (e->stepnum < stepnum) {
1561       ierr = TurnForward(ts);CHKERRQ(ierr);
1562       ierr = ReCompute(ts,tjsch,e->stepnum,stepnum);CHKERRQ(ierr);
1563     }
1564     if (e->stepnum == stepnum) {
1565       ierr = StackPop(stack,&e);CHKERRQ(ierr);
1566     }
1567   }
1568   tjsch->rctx->reverseonestep = PETSC_FALSE;
1569   PetscFunctionReturn(0);
1570 }
1571 
1572 static PetscErrorCode SetTrajRMS(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X)
1573 {
1574   Stack           *stack = &tjsch->stack;
1575   PetscInt        store;
1576   StackElement    e;
1577   PetscRevolveInt rtotal_steps,rstepnum;
1578   PetscErrorCode  ierr;
1579 
1580   PetscFunctionBegin;
1581   if (!stack->solution_only && stepnum == 0) PetscFunctionReturn(0);
1582   if (stack->solution_only && stepnum == tjsch->total_steps) PetscFunctionReturn(0);
1583   ierr = PetscRevolveIntCast(tjsch->total_steps,&rtotal_steps);CHKERRQ(ierr);
1584   ierr = PetscRevolveIntCast(stepnum,&rstepnum);CHKERRQ(ierr);
1585   ierr = ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,rtotal_steps,rstepnum,rstepnum,PETSC_FALSE,&store);CHKERRQ(ierr);
1586   if (store == 1) {
1587     CheckpointType cptype;
1588     if (stepnum < stack->top) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_MEMC,"Illegal modification of a non-top stack element");
1589     cptype = stack->solution_only ? SOLUTIONONLY : SOLUTION_STAGES;
1590     ierr = ElementCreate(ts,cptype,stack,&e);CHKERRQ(ierr);
1591     ierr = ElementSet(ts,stack,&e,stepnum,time,X);CHKERRQ(ierr);
1592     ierr = StackPush(stack,e);CHKERRQ(ierr);
1593   } else if (store == 2) {
1594     ierr = DumpSingle(tj,ts,stack,tjsch->rctx->check+1);CHKERRQ(ierr);
1595   }
1596   PetscFunctionReturn(0);
1597 }
1598 
1599 static PetscErrorCode GetTrajRMS(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum)
1600 {
1601   Stack           *stack = &tjsch->stack;
1602   PetscRevolveInt whattodo,shift,rstepnum;
1603   PetscInt        restart;
1604   PetscBool       ondisk;
1605   StackElement    e;
1606   PetscErrorCode  ierr;
1607 
1608   PetscFunctionBegin;
1609   if (stepnum == 0 || stepnum == tjsch->total_steps) {
1610     ierr = TurnBackward(ts);CHKERRQ(ierr);
1611     tjsch->rctx->reverseonestep = PETSC_FALSE;
1612     PetscFunctionReturn(0);
1613   }
1614   ierr = PetscRevolveIntCast(stepnum,&rstepnum);CHKERRQ(ierr);
1615   tjsch->rctx->capo = rstepnum;
1616   tjsch->rctx->oldcapo = tjsch->rctx->capo;
1617   shift = 0;
1618   whattodo = revolve_action(&tjsch->rctx->check,&tjsch->rctx->capo,&tjsch->rctx->fine,tjsch->rctx->snaps_in,&tjsch->rctx->info,&tjsch->rctx->where); /* whattodo=restore */
1619   ierr = printwhattodo(tj->monitor,whattodo,tjsch->rctx,shift);CHKERRQ(ierr);
1620   /* restore a checkpoint */
1621   restart = tjsch->rctx->capo;
1622   if (!tjsch->rctx->where) {
1623     ondisk = PETSC_TRUE;
1624     ierr = LoadSingle(tj,ts,stack,tjsch->rctx->check+1);CHKERRQ(ierr);
1625     ierr = TurnBackward(ts);CHKERRQ(ierr);
1626   } else {
1627     ondisk = PETSC_FALSE;
1628     ierr = StackTop(stack,&e);CHKERRQ(ierr);
1629     ierr = UpdateTS(ts,stack,e,stepnum,PETSC_TRUE);CHKERRQ(ierr);
1630   }
1631   if (!stack->solution_only) { /* whattodo must be 5 or 8 */
1632     /* ask Revolve what to do next */
1633     tjsch->rctx->oldcapo = tjsch->rctx->capo;
1634     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*/
1635     ierr = printwhattodo(tj->monitor,whattodo,tjsch->rctx,shift);CHKERRQ(ierr);
1636     if (whattodo == 3 || whattodo == 4) tjsch->rctx->reverseonestep = PETSC_TRUE;
1637     if (whattodo == 1) tjsch->rctx->stepsleft = tjsch->rctx->capo-tjsch->rctx->oldcapo;
1638     if (tj->monitor) {
1639       ierr = PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr);
1640       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);
1641       ierr = PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr);
1642     }
1643     if (!tjsch->rctx->reverseonestep && tjsch->rctx->stepsleft > 0) tjsch->rctx->stepsleft--;
1644     restart++; /* skip one step */
1645   }
1646   if (stack->solution_only || (!stack->solution_only && restart < stepnum)) {
1647     ierr = TurnForward(ts);CHKERRQ(ierr);
1648     ierr = ReCompute(ts,tjsch,restart,stepnum);CHKERRQ(ierr);
1649   }
1650   if (!ondisk && ( (stack->solution_only && e->stepnum+1 == stepnum) || (!stack->solution_only && e->stepnum == stepnum))) {
1651     ierr = StackPop(stack,&e);CHKERRQ(ierr);
1652   }
1653   tjsch->rctx->reverseonestep = PETSC_FALSE;
1654   PetscFunctionReturn(0);
1655 }
1656 #endif
1657 
1658 #if defined(PETSC_HAVE_CAMS)
1659 /* Optimal offline adjoint checkpointing for multistage time integration methods */
1660 static PetscErrorCode SetTrajAOF(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X)
1661 {
1662   Stack          *stack = &tjsch->stack;
1663   StackElement   e;
1664   PetscErrorCode ierr;
1665 
1666   PetscFunctionBegin;
1667   /* skip if no checkpoint to use. This also avoids an error when num_units_avail=0  */
1668   if (tjsch->actx->nextcheckpointstep == -1) PetscFunctionReturn(0);
1669   if (stepnum == 0) { /* When placing the first checkpoint, no need to change the units available */
1670     if (stack->solution_only) {
1671       ierr = offline_ca(tjsch->actx->lastcheckpointstep,tjsch->actx->num_units_avail,tjsch->actx->endstep,&tjsch->actx->nextcheckpointstep);CHKERRQ(ierr);
1672     } else {
1673       /* First two arguments must be -1 when first time calling cams */
1674       ierr = offline_cams(tjsch->actx->lastcheckpointstep,tjsch->actx->lastcheckpointtype,tjsch->actx->num_units_avail,tjsch->actx->endstep,tjsch->actx->num_stages,&tjsch->actx->nextcheckpointstep,&tjsch->actx->nextcheckpointtype);CHKERRQ(ierr);
1675     }
1676   }
1677 
1678   if (stack->solution_only && stepnum == tjsch->total_steps) PetscFunctionReturn(0);
1679 
1680   if (tjsch->actx->nextcheckpointstep == stepnum) {
1681     if (stepnum < stack->top) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_MEMC,"Illegal modification of a non-top stack element");
1682 
1683     if (tjsch->actx->nextcheckpointtype == 2) { /* solution + stage values */
1684       if (tj->monitor) {
1685         ierr = PetscViewerASCIIPrintf(tj->monitor,"Store in checkpoint number %D with stage values and solution (located in RAM)\n",stepnum);CHKERRQ(ierr);
1686       }
1687       ierr = ElementCreate(ts,SOLUTION_STAGES,stack,&e);CHKERRQ(ierr);
1688       ierr = ElementSet(ts,stack,&e,stepnum,time,X);CHKERRQ(ierr);
1689     }
1690     if (tjsch->actx->nextcheckpointtype == 1) {
1691       if (tj->monitor) {
1692         ierr = PetscViewerASCIIPrintf(tj->monitor,"Store in checkpoint number %D with stage values (located in RAM)\n",stepnum);CHKERRQ(ierr);
1693       }
1694       ierr = ElementCreate(ts,STAGESONLY,stack,&e);CHKERRQ(ierr);
1695       ierr = ElementSet(ts,stack,&e,stepnum,time,X);CHKERRQ(ierr);
1696     }
1697     if (tjsch->actx->nextcheckpointtype == 0) { /* solution only */
1698       if (tj->monitor) {
1699         ierr = PetscViewerASCIIPrintf(tj->monitor,"Store in checkpoint number %D (located in RAM)\n",stepnum);CHKERRQ(ierr);
1700       }
1701       ierr = ElementCreate(ts,SOLUTIONONLY,stack,&e);CHKERRQ(ierr);
1702       ierr = ElementSet(ts,stack,&e,stepnum,time,X);CHKERRQ(ierr);
1703     }
1704     ierr = StackPush(stack,e);CHKERRQ(ierr);
1705 
1706     tjsch->actx->lastcheckpointstep = stepnum;
1707     if (stack->solution_only) {
1708       ierr = offline_ca(tjsch->actx->lastcheckpointstep,tjsch->actx->num_units_avail,tjsch->actx->endstep,&tjsch->actx->nextcheckpointstep);CHKERRQ(ierr);
1709       tjsch->actx->num_units_avail--;
1710     } else {
1711       tjsch->actx->lastcheckpointtype = tjsch->actx->nextcheckpointtype;
1712       ierr = offline_cams(tjsch->actx->lastcheckpointstep,tjsch->actx->lastcheckpointtype,tjsch->actx->num_units_avail,tjsch->actx->endstep,tjsch->actx->num_stages,&tjsch->actx->nextcheckpointstep,&tjsch->actx->nextcheckpointtype);CHKERRQ(ierr);
1713       if (tjsch->actx->lastcheckpointtype == 2) tjsch->actx->num_units_avail -= tjsch->actx->num_stages+1;
1714       if (tjsch->actx->lastcheckpointtype == 1) tjsch->actx->num_units_avail -= tjsch->actx->num_stages;
1715       if (tjsch->actx->lastcheckpointtype == 0) tjsch->actx->num_units_avail--;
1716     }
1717   }
1718   PetscFunctionReturn(0);
1719 }
1720 
1721 static PetscErrorCode GetTrajAOF(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum)
1722 {
1723   Stack          *stack = &tjsch->stack;
1724   StackElement   e;
1725   PetscInt       estepnum;
1726   PetscErrorCode ierr;
1727 
1728   PetscFunctionBegin;
1729   if (stepnum == 0 || stepnum == tjsch->total_steps) {
1730     ierr = TurnBackward(ts);CHKERRQ(ierr);
1731     PetscFunctionReturn(0);
1732   }
1733   /* Restore a checkpoint */
1734   ierr = StackTop(stack,&e);CHKERRQ(ierr);
1735   estepnum = e->stepnum;
1736   if (estepnum == stepnum && e->cptype == SOLUTIONONLY) { /* discard the checkpoint if not useful (corner case) */
1737     ierr = StackPop(stack,&e);CHKERRQ(ierr);
1738     tjsch->actx->num_units_avail++;
1739     ierr = StackTop(stack,&e);CHKERRQ(ierr);
1740     estepnum = e->stepnum;
1741   }
1742   /* Update TS with stage values if an adjoint step can be taken immediately */
1743   if (HaveStages(e->cptype)) {
1744     if (tj->monitor) {
1745       ierr = PetscViewerASCIIPrintf(tj->monitor,"Restore in checkpoint number %D with stage values (located in RAM)\n",e->stepnum);CHKERRQ(ierr);
1746     }
1747     if (e->cptype == STAGESONLY) tjsch->actx->num_units_avail += tjsch->actx->num_stages;
1748     if (e->cptype == SOLUTION_STAGES) tjsch->actx->num_units_avail += tjsch->actx->num_stages+1;
1749   } else {
1750     if (tj->monitor) {
1751       ierr = PetscViewerASCIIPrintf(tj->monitor,"Restore in checkpoint number %D (located in RAM)\n",e->stepnum);CHKERRQ(ierr);
1752     }
1753     tjsch->actx->num_units_avail++;
1754   }
1755   ierr = UpdateTS(ts,stack,e,stepnum,PETSC_TRUE);CHKERRQ(ierr);
1756   /* Query the scheduler */
1757   tjsch->actx->lastcheckpointstep = estepnum;
1758   tjsch->actx->endstep = stepnum;
1759   if (stack->solution_only) { /* start with restoring a checkpoint */
1760     ierr = offline_ca(tjsch->actx->lastcheckpointstep,tjsch->actx->num_units_avail,tjsch->actx->endstep,&tjsch->actx->nextcheckpointstep);CHKERRQ(ierr);
1761   } else { /* 2 revolve actions: restore a checkpoint and then advance */
1762     tjsch->actx->lastcheckpointtype = e->cptype;
1763     ierr = offline_cams(tjsch->actx->lastcheckpointstep,tjsch->actx->lastcheckpointtype,tjsch->actx->num_units_avail,tjsch->actx->endstep,tjsch->actx->num_stages,&tjsch->actx->nextcheckpointstep, &tjsch->actx->nextcheckpointtype);CHKERRQ(ierr);
1764   }
1765   /* Discard the checkpoint if not needed, decrease the number of available checkpoints if it still stays in stack */
1766   if (HaveStages(e->cptype)) {
1767     if (estepnum == stepnum) {
1768       ierr = StackPop(stack,&e);CHKERRQ(ierr);
1769     } else {
1770       if (e->cptype == STAGESONLY) tjsch->actx->num_units_avail -= tjsch->actx->num_stages;
1771       if (e->cptype == SOLUTION_STAGES) tjsch->actx->num_units_avail -= tjsch->actx->num_stages+1;
1772     }
1773   } else {
1774     if (estepnum+1 == stepnum) {
1775       ierr = StackPop(stack,&e);CHKERRQ(ierr);
1776     } else {
1777       tjsch->actx->num_units_avail--;
1778     }
1779   }
1780   /* Recompute from the restored checkpoint */
1781   if (stack->solution_only || (!stack->solution_only && estepnum < stepnum)) {
1782     ierr = TurnForward(ts);CHKERRQ(ierr);
1783     ierr = ReCompute(ts,tjsch,estepnum,stepnum);CHKERRQ(ierr);
1784   }
1785   PetscFunctionReturn(0);
1786 }
1787 #endif
1788 
1789 static PetscErrorCode TSTrajectorySet_Memory(TSTrajectory tj,TS ts,PetscInt stepnum,PetscReal time,Vec X)
1790 {
1791   TJScheduler *tjsch = (TJScheduler*)tj->data;
1792   PetscErrorCode ierr;
1793 
1794   PetscFunctionBegin;
1795   if (!tjsch->recompute) { /* use global stepnum in the forward sweep */
1796     ierr = TSGetStepNumber(ts,&stepnum);CHKERRQ(ierr);
1797   }
1798   /* for consistency */
1799   if (!tjsch->recompute && stepnum == 0) ts->ptime_prev = ts->ptime-ts->time_step;
1800   switch (tjsch->stype) {
1801     case NONE:
1802       if (tj->adjoint_solve_mode) {
1803         ierr = SetTrajN(ts,tjsch,stepnum,time,X);CHKERRQ(ierr);
1804       } else {
1805         ierr = SetTrajN_2(ts,tjsch,stepnum,time,X);CHKERRQ(ierr);
1806       }
1807       break;
1808     case TWO_LEVEL_NOREVOLVE:
1809       if (!tj->adjoint_solve_mode) SETERRQ(PetscObjectComm((PetscObject)tj),PETSC_ERR_SUP,"Not implemented");
1810       ierr = SetTrajTLNR(tj,ts,tjsch,stepnum,time,X);CHKERRQ(ierr);
1811       break;
1812 #if defined(PETSC_HAVE_REVOLVE)
1813     case TWO_LEVEL_REVOLVE:
1814       if (!tj->adjoint_solve_mode) SETERRQ(PetscObjectComm((PetscObject)tj),PETSC_ERR_SUP,"Not implemented");
1815       ierr = SetTrajTLR(tj,ts,tjsch,stepnum,time,X);CHKERRQ(ierr);
1816       break;
1817     case TWO_LEVEL_TWO_REVOLVE:
1818       if (!tj->adjoint_solve_mode) SETERRQ(PetscObjectComm((PetscObject)tj),PETSC_ERR_SUP,"Not implemented");
1819       ierr = SetTrajTLTR(tj,ts,tjsch,stepnum,time,X);CHKERRQ(ierr);
1820       break;
1821     case REVOLVE_OFFLINE:
1822       if (!tj->adjoint_solve_mode) SETERRQ(PetscObjectComm((PetscObject)tj),PETSC_ERR_SUP,"Not implemented");
1823       ierr = SetTrajROF(tj,ts,tjsch,stepnum,time,X);CHKERRQ(ierr);
1824       break;
1825     case REVOLVE_ONLINE:
1826       if (!tj->adjoint_solve_mode) SETERRQ(PetscObjectComm((PetscObject)tj),PETSC_ERR_SUP,"Not implemented");
1827       ierr = SetTrajRON(tj,ts,tjsch,stepnum,time,X);CHKERRQ(ierr);
1828       break;
1829     case REVOLVE_MULTISTAGE:
1830       if (!tj->adjoint_solve_mode) SETERRQ(PetscObjectComm((PetscObject)tj),PETSC_ERR_SUP,"Not implemented");
1831       ierr = SetTrajRMS(tj,ts,tjsch,stepnum,time,X);CHKERRQ(ierr);
1832       break;
1833 #endif
1834 #if defined(PETSC_HAVE_CAMS)
1835     case CAMS_OFFLINE:
1836       if (!tj->adjoint_solve_mode) SETERRQ(PetscObjectComm((PetscObject)tj),PETSC_ERR_SUP,"Not implemented");
1837       ierr = SetTrajAOF(tj,ts,tjsch,stepnum,time,X);CHKERRQ(ierr);
1838       break;
1839 #endif
1840     default:
1841       break;
1842   }
1843   PetscFunctionReturn(0);
1844 }
1845 
1846 static PetscErrorCode TSTrajectoryGet_Memory(TSTrajectory tj,TS ts,PetscInt stepnum,PetscReal *t)
1847 {
1848   TJScheduler *tjsch = (TJScheduler*)tj->data;
1849   PetscErrorCode ierr;
1850 
1851   PetscFunctionBegin;
1852   if (tj->adjoint_solve_mode && stepnum == 0) {
1853     ierr = TSTrajectoryReset(tj);CHKERRQ(ierr); /* reset TSTrajectory so users do not need to reset TSTrajectory */
1854     PetscFunctionReturn(0);
1855   }
1856   switch (tjsch->stype) {
1857     case NONE:
1858       if (tj->adjoint_solve_mode) {
1859         ierr = GetTrajN(ts,tjsch,stepnum);CHKERRQ(ierr);
1860       } else {
1861         ierr = GetTrajN_2(ts,tjsch,stepnum);CHKERRQ(ierr);
1862       }
1863       break;
1864     case TWO_LEVEL_NOREVOLVE:
1865       if (!tj->adjoint_solve_mode) SETERRQ(PetscObjectComm((PetscObject)tj),PETSC_ERR_SUP,"Not implemented");
1866       ierr = GetTrajTLNR(tj,ts,tjsch,stepnum);CHKERRQ(ierr);
1867       break;
1868 #if defined(PETSC_HAVE_REVOLVE)
1869     case TWO_LEVEL_REVOLVE:
1870       if (!tj->adjoint_solve_mode) SETERRQ(PetscObjectComm((PetscObject)tj),PETSC_ERR_SUP,"Not implemented");
1871       ierr = GetTrajTLR(tj,ts,tjsch,stepnum);CHKERRQ(ierr);
1872       break;
1873     case TWO_LEVEL_TWO_REVOLVE:
1874       if (!tj->adjoint_solve_mode) SETERRQ(PetscObjectComm((PetscObject)tj),PETSC_ERR_SUP,"Not implemented");
1875       ierr = GetTrajTLTR(tj,ts,tjsch,stepnum);CHKERRQ(ierr);
1876       break;
1877     case REVOLVE_OFFLINE:
1878       if (!tj->adjoint_solve_mode) SETERRQ(PetscObjectComm((PetscObject)tj),PETSC_ERR_SUP,"Not implemented");
1879       ierr = GetTrajROF(tj,ts,tjsch,stepnum);CHKERRQ(ierr);
1880       break;
1881     case REVOLVE_ONLINE:
1882       if (!tj->adjoint_solve_mode) SETERRQ(PetscObjectComm((PetscObject)tj),PETSC_ERR_SUP,"Not implemented");
1883       ierr = GetTrajRON(tj,ts,tjsch,stepnum);CHKERRQ(ierr);
1884       break;
1885     case REVOLVE_MULTISTAGE:
1886       if (!tj->adjoint_solve_mode) SETERRQ(PetscObjectComm((PetscObject)tj),PETSC_ERR_SUP,"Not implemented");
1887       ierr = GetTrajRMS(tj,ts,tjsch,stepnum);CHKERRQ(ierr);
1888       break;
1889 #endif
1890 #if defined(PETSC_HAVE_CAMS)
1891     case CAMS_OFFLINE:
1892       if (!tj->adjoint_solve_mode) SETERRQ(PetscObjectComm((PetscObject)tj),PETSC_ERR_SUP,"Not implemented");
1893       ierr = GetTrajAOF(tj,ts,tjsch,stepnum);CHKERRQ(ierr);
1894       break;
1895 #endif
1896     default:
1897       break;
1898   }
1899   PetscFunctionReturn(0);
1900 }
1901 
1902 PETSC_UNUSED static PetscErrorCode TSTrajectorySetStride_Memory(TSTrajectory tj,PetscInt stride)
1903 {
1904   TJScheduler *tjsch = (TJScheduler*)tj->data;
1905 
1906   PetscFunctionBegin;
1907   tjsch->stride = stride;
1908   PetscFunctionReturn(0);
1909 }
1910 
1911 static PetscErrorCode TSTrajectorySetMaxCpsRAM_Memory(TSTrajectory tj,PetscInt max_cps_ram)
1912 {
1913   TJScheduler *tjsch = (TJScheduler*)tj->data;
1914 
1915   PetscFunctionBegin;
1916   tjsch->max_cps_ram = max_cps_ram;
1917   PetscFunctionReturn(0);
1918 }
1919 
1920 static PetscErrorCode TSTrajectorySetMaxCpsDisk_Memory(TSTrajectory tj,PetscInt max_cps_disk)
1921 {
1922   TJScheduler *tjsch = (TJScheduler*)tj->data;
1923 
1924   PetscFunctionBegin;
1925   tjsch->max_cps_disk = max_cps_disk;
1926   PetscFunctionReturn(0);
1927 }
1928 
1929 static PetscErrorCode TSTrajectorySetMaxUnitsRAM_Memory(TSTrajectory tj,PetscInt max_units_ram)
1930 {
1931   TJScheduler *tjsch = (TJScheduler*)tj->data;
1932 
1933   PetscFunctionBegin;
1934   if (!tjsch->max_cps_ram) SETERRQ(PetscObjectComm((PetscObject)tj),PETSC_ERR_ARG_INCOMP,"Conflict with -ts_trjaectory_max_cps_ram or TSTrajectorySetMaxCpsRAM. You can set max_cps_ram or max_units_ram, but not both at the same time.");
1935   tjsch->max_units_ram = max_units_ram;
1936   PetscFunctionReturn(0);
1937 }
1938 
1939 static PetscErrorCode TSTrajectorySetMaxUnitsDisk_Memory(TSTrajectory tj,PetscInt max_units_disk)
1940 {
1941   TJScheduler *tjsch = (TJScheduler*)tj->data;
1942 
1943   PetscFunctionBegin;
1944   if (!tjsch->max_cps_disk) SETERRQ(PetscObjectComm((PetscObject)tj),PETSC_ERR_ARG_INCOMP,"Conflict with -ts_trjaectory_max_cps_disk or TSTrajectorySetMaxCpsDisk. You can set max_cps_disk or max_units_disk, but not both at the same time.");
1945   tjsch->max_units_ram = max_units_disk;
1946   PetscFunctionReturn(0);
1947 }
1948 
1949 #if defined(PETSC_HAVE_REVOLVE)
1950 PETSC_UNUSED static PetscErrorCode TSTrajectorySetRevolveOnline(TSTrajectory tj,PetscBool use_online)
1951 {
1952   TJScheduler *tjsch = (TJScheduler*)tj->data;
1953 
1954   PetscFunctionBegin;
1955   tjsch->use_online = use_online;
1956   PetscFunctionReturn(0);
1957 }
1958 #endif
1959 
1960 PETSC_UNUSED static PetscErrorCode TSTrajectorySetSaveStack(TSTrajectory tj,PetscBool save_stack)
1961 {
1962   TJScheduler *tjsch = (TJScheduler*)tj->data;
1963 
1964   PetscFunctionBegin;
1965   tjsch->save_stack = save_stack;
1966   PetscFunctionReturn(0);
1967 }
1968 
1969 PETSC_UNUSED static PetscErrorCode TSTrajectorySetUseDRAM(TSTrajectory tj,PetscBool use_dram)
1970 {
1971   TJScheduler *tjsch = (TJScheduler*)tj->data;
1972 
1973   PetscFunctionBegin;
1974   tjsch->stack.use_dram = use_dram;
1975   PetscFunctionReturn(0);
1976 }
1977 
1978 /*@C
1979   TSTrajectorySetMaxCpsRAM - Set maximum number of checkpoints in RAM
1980 
1981   Logically collective
1982 
1983   Input Parameter:
1984 .  tj - tstrajectory context
1985 
1986   Output Parameter:
1987 .  max_cps_ram - maximum number of checkpoints in RAM
1988 
1989   Level: intermediate
1990 
1991 .seealso: TSTrajectorySetMaxUnitsRAM()
1992 @*/
1993 PetscErrorCode TSTrajectorySetMaxCpsRAM(TSTrajectory tj,PetscInt max_cps_ram)
1994 {
1995   PetscErrorCode ierr;
1996 
1997   PetscFunctionBegin;
1998   ierr = PetscUseMethod(tj,"TSTrajectorySetMaxCpsRAM_C",(TSTrajectory,PetscInt),(tj,max_cps_ram));CHKERRQ(ierr);
1999   PetscFunctionReturn(0);
2000 }
2001 
2002 /*@C
2003   TSTrajectorySetMaxCpsDisk - Set maximum number of checkpoints on disk
2004 
2005   Logically collective
2006 
2007   Input Parameter:
2008 .  tj - tstrajectory context
2009 
2010   Output Parameter:
2011 .  max_cps_disk - maximum number of checkpoints on disk
2012 
2013   Level: intermediate
2014 
2015 .seealso: TSTrajectorySetMaxUnitsDisk(), TSTrajectorySetMaxUnitsRAM()
2016 @*/
2017 PetscErrorCode TSTrajectorySetMaxCpsDisk(TSTrajectory tj,PetscInt max_cps_disk)
2018 {
2019   PetscErrorCode ierr;
2020 
2021   PetscFunctionBegin;
2022   ierr = PetscUseMethod(tj,"TSTrajectorySetMaxCpsDisk_C",(TSTrajectory,PetscInt),(tj,max_cps_disk));CHKERRQ(ierr);
2023   PetscFunctionReturn(0);
2024 }
2025 
2026 /*@C
2027   TSTrajectorySetMaxUnitsRAM - Set maximum number of checkpointing units in RAM
2028 
2029   Logically collective
2030 
2031   Input Parameter:
2032 .  tj - tstrajectory context
2033 
2034   Output Parameter:
2035 .  max_units_ram - maximum number of checkpointing units in RAM
2036 
2037   Level: intermediate
2038 
2039 .seealso: TSTrajectorySetMaxCpsRAM()
2040 @*/
2041 PetscErrorCode TSTrajectorySetMaxUnitsRAM(TSTrajectory tj,PetscInt max_units_ram)
2042 {
2043   PetscErrorCode ierr;
2044 
2045   PetscFunctionBegin;
2046   ierr = PetscUseMethod(tj,"TSTrajectorySetMaxUnitsRAM_C",(TSTrajectory,PetscInt),(tj,max_units_ram));CHKERRQ(ierr);
2047   PetscFunctionReturn(0);
2048 }
2049 
2050 /*@C
2051   TSTrajectorySetMaxUnitsDisk - Set maximum number of checkpointing units on disk
2052 
2053   Logically collective
2054 
2055   Input Parameter:
2056 .  tj - tstrajectory context
2057 
2058   Output Parameter:
2059 .  max_units_disk - maximum number of checkpointing units on disk
2060 
2061   Level: intermediate
2062 
2063 .seealso: TSTrajectorySetMaxCpsDisk()
2064 @*/
2065 PetscErrorCode TSTrajectorySetMaxUnitsDisk(TSTrajectory tj,PetscInt max_units_disk)
2066 {
2067   PetscErrorCode ierr;
2068 
2069   PetscFunctionBegin;
2070   ierr = PetscUseMethod(tj,"TSTrajectorySetMaxUnitsDisk_C",(TSTrajectory,PetscInt),(tj,max_units_disk));CHKERRQ(ierr);
2071   PetscFunctionReturn(0);
2072 }
2073 
2074 static PetscErrorCode TSTrajectorySetFromOptions_Memory(PetscOptionItems *PetscOptionsObject,TSTrajectory tj)
2075 {
2076   TJScheduler    *tjsch = (TJScheduler*)tj->data;
2077   PetscInt       max_cps_ram,max_cps_disk,max_units_ram,max_units_disk;
2078   PetscBool      flg;
2079   PetscErrorCode ierr;
2080 
2081   PetscFunctionBegin;
2082   ierr = PetscOptionsHead(PetscOptionsObject,"Memory based TS trajectory options");CHKERRQ(ierr);
2083   {
2084     ierr = PetscOptionsInt("-ts_trajectory_max_cps_ram","Maximum number of checkpoints in RAM","TSTrajectorySetMaxCpsRAM",tjsch->max_cps_ram,&max_cps_ram,&flg);CHKERRQ(ierr);
2085     if (flg) {
2086       ierr = TSTrajectorySetMaxCpsRAM(tj,max_cps_ram);CHKERRQ(ierr);
2087     }
2088     ierr = PetscOptionsInt("-ts_trajectory_max_cps_disk","Maximum number of checkpoints on disk","TSTrajectorySetMaxCpsDisk",tjsch->max_cps_disk,&max_cps_disk,&flg);CHKERRQ(ierr);
2089     if (flg) {
2090       ierr = TSTrajectorySetMaxCpsDisk(tj,max_cps_disk);CHKERRQ(ierr);
2091     }
2092     ierr = PetscOptionsInt("-ts_trajectory_max_units_ram","Maximum number of checkpointing units in RAM","TSTrajectorySetMaxUnitsRAM",tjsch->max_units_ram,&max_units_ram,&flg);CHKERRQ(ierr);
2093     if (flg) {
2094       ierr = TSTrajectorySetMaxUnitsRAM(tj,max_units_ram);CHKERRQ(ierr);
2095     }
2096     ierr = PetscOptionsInt("-ts_trajectory_max_units_disk","Maximum number of checkpointing units on disk","TSTrajectorySetMaxUnitsDisk",tjsch->max_units_disk,&max_units_disk,&flg);CHKERRQ(ierr);
2097     if (flg) {
2098       ierr = TSTrajectorySetMaxUnitsDisk(tj,max_units_disk);CHKERRQ(ierr);
2099     }
2100     ierr = PetscOptionsInt("-ts_trajectory_stride","Stride to save checkpoints to file","TSTrajectorySetStride",tjsch->stride,&tjsch->stride,NULL);CHKERRQ(ierr);
2101 #if defined(PETSC_HAVE_REVOLVE)
2102     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);
2103 #endif
2104     ierr = PetscOptionsBool("-ts_trajectory_save_stack","Save all stack to disk","TSTrajectorySetSaveStack",tjsch->save_stack,&tjsch->save_stack,NULL);CHKERRQ(ierr);
2105     ierr = PetscOptionsBool("-ts_trajectory_use_dram","Use DRAM for checkpointing","TSTrajectorySetUseDRAM",tjsch->stack.use_dram,&tjsch->stack.use_dram,NULL);CHKERRQ(ierr);
2106   }
2107   ierr = PetscOptionsTail();CHKERRQ(ierr);
2108   tjsch->stack.solution_only = tj->solution_only;
2109   PetscFunctionReturn(0);
2110 }
2111 
2112 static PetscErrorCode TSTrajectorySetUp_Memory(TSTrajectory tj,TS ts)
2113 {
2114   TJScheduler    *tjsch = (TJScheduler*)tj->data;
2115   Stack          *stack = &tjsch->stack;
2116 #if defined(PETSC_HAVE_REVOLVE)
2117   RevolveCTX     *rctx,*rctx2;
2118   DiskStack      *diskstack = &tjsch->diskstack;
2119   PetscInt       diskblocks;
2120 #endif
2121   PetscInt       numY,total_steps;
2122   PetscBool      fixedtimestep;
2123   PetscErrorCode ierr;
2124 
2125   PetscFunctionBegin;
2126   if (ts->adapt) {
2127     ierr = PetscObjectTypeCompare((PetscObject)ts->adapt,TSADAPTNONE,&fixedtimestep);CHKERRQ(ierr);
2128   } else {
2129     fixedtimestep = PETSC_TRUE;
2130   }
2131   total_steps = (PetscInt)(PetscCeilReal((ts->max_time-ts->ptime)/ts->time_step));
2132   total_steps = total_steps < 0 ? PETSC_MAX_INT : total_steps;
2133   if (fixedtimestep) tjsch->total_steps = PetscMin(ts->max_steps,total_steps);
2134 
2135   ierr = TSGetStages(ts,&numY,PETSC_IGNORE);CHKERRQ(ierr);
2136   if (stack->solution_only) {
2137     if (tjsch->max_units_ram) tjsch->max_cps_ram = tjsch->max_units_ram;
2138     else tjsch->max_units_ram = tjsch->max_cps_ram;
2139     if (tjsch->max_units_disk) tjsch->max_cps_disk = tjsch->max_units_disk;
2140   } else {
2141     if (tjsch->max_units_ram) tjsch->max_cps_ram = (ts->stifflyaccurate) ? tjsch->max_units_ram/numY : tjsch->max_units_ram/(numY+1);
2142     else tjsch->max_units_ram = (ts->stifflyaccurate) ? numY*tjsch->max_cps_ram : (numY+1)*tjsch->max_cps_ram;
2143     if (tjsch->max_units_disk) tjsch->max_cps_disk = (ts->stifflyaccurate) ? tjsch->max_units_disk/numY : tjsch->max_units_disk/(numY+1);
2144     else tjsch->max_units_disk = (ts->stifflyaccurate) ? numY*tjsch->max_cps_disk : (numY+1)*tjsch->max_cps_disk;
2145   }
2146   if (tjsch->max_cps_ram > 0) stack->stacksize = tjsch->max_units_ram; /* maximum stack size. Could be overallocated. */
2147 
2148   /* Determine the scheduler type */
2149   if (tjsch->stride > 1) { /* two level mode */
2150     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.");
2151     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 */
2152     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 */
2153     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 */
2154   } else { /* single level mode */
2155     if (fixedtimestep) {
2156       if (tjsch->max_cps_ram >= tjsch->total_steps-1 || tjsch->max_cps_ram == -1)
2157         tjsch->stype = NONE; /* checkpoint all */
2158       else
2159 #if defined(PETSC_HAVE_CAMS)
2160       /* CAMS has the priority over revolve */
2161         tjsch->stype = CAMS_OFFLINE;
2162 #else
2163         tjsch->stype = (tjsch->max_cps_disk>1) ? REVOLVE_MULTISTAGE : REVOLVE_OFFLINE;
2164 #endif
2165     } else tjsch->stype = NONE; /* checkpoint all for adaptive time step */
2166 #if defined(PETSC_HAVE_REVOLVE)
2167     if (tjsch->use_online) tjsch->stype = REVOLVE_ONLINE; /* trick into online (for testing purpose only) */
2168 #endif
2169     if (tjsch->stype != NONE && tjsch->max_cps_ram < 1 && tjsch->max_cps_disk < 1) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_INCOMP,"The specified storage capacity is insufficient for one checkpoint, which is the minimum");
2170   }
2171   if (tjsch->stype >= CAMS_OFFLINE) {
2172 #ifndef PETSC_HAVE_CAMS
2173     SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"CAMS 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-cams.");
2174 #else
2175     CAMSCTX  *actx;
2176     PetscInt ns = 0;
2177     if (stack->solution_only) {
2178       offline_ca_create(tjsch->total_steps,tjsch->max_cps_ram);
2179     } else {
2180       ierr = TSGetStages(ts,&ns,PETSC_IGNORE);CHKERRQ(ierr);
2181       offline_cams_create(tjsch->total_steps,tjsch->max_units_ram,ns,ts->stifflyaccurate);
2182     }
2183     ierr = PetscNew(&actx);CHKERRQ(ierr);
2184     actx->lastcheckpointstep    = -1; /* -1 can trigger the initialization of CAMS */
2185     actx->lastcheckpointtype    = -1; /* -1 can trigger the initialization of CAMS */
2186     actx->endstep               = tjsch->total_steps;
2187     actx->num_units_avail       = tjsch->max_units_ram;
2188     actx->num_stages            = ns;
2189     tjsch->actx                 = actx;
2190 #endif
2191   } else if (tjsch->stype > TWO_LEVEL_NOREVOLVE) {
2192 #ifndef PETSC_HAVE_REVOLVE
2193     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.");
2194 #else
2195     PetscRevolveInt rfine,rsnaps,rsnaps2;
2196 
2197     switch (tjsch->stype) {
2198       case TWO_LEVEL_REVOLVE:
2199         ierr = PetscRevolveIntCast(tjsch->stride,&rfine);CHKERRQ(ierr);
2200         ierr = PetscRevolveIntCast(tjsch->max_cps_ram,&rsnaps);CHKERRQ(ierr);
2201         revolve_create_offline(rfine,rsnaps);
2202         break;
2203       case TWO_LEVEL_TWO_REVOLVE:
2204         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. */
2205         diskstack->stacksize = diskblocks;
2206         ierr = PetscRevolveIntCast(tjsch->stride,&rfine);CHKERRQ(ierr);
2207         ierr = PetscRevolveIntCast(tjsch->max_cps_ram,&rsnaps);CHKERRQ(ierr);
2208         revolve_create_offline(rfine,rsnaps);
2209         ierr = PetscRevolveIntCast((tjsch->total_steps+tjsch->stride-1)/tjsch->stride,&rfine);CHKERRQ(ierr);
2210         ierr = PetscRevolveIntCast(diskblocks,&rsnaps);CHKERRQ(ierr);
2211         revolve2_create_offline(rfine,rsnaps);
2212         ierr = PetscNew(&rctx2);CHKERRQ(ierr);
2213         rctx2->snaps_in       = rsnaps;
2214         rctx2->reverseonestep = PETSC_FALSE;
2215         rctx2->check          = 0;
2216         rctx2->oldcapo        = 0;
2217         rctx2->capo           = 0;
2218         rctx2->info           = 2;
2219         rctx2->fine           = rfine;
2220         tjsch->rctx2          = rctx2;
2221         diskstack->top        = -1;
2222         ierr = PetscMalloc1(diskstack->stacksize*sizeof(PetscInt),&diskstack->container);CHKERRQ(ierr);
2223         break;
2224       case REVOLVE_OFFLINE:
2225         ierr = PetscRevolveIntCast(tjsch->total_steps,&rfine);CHKERRQ(ierr);
2226         ierr = PetscRevolveIntCast(tjsch->max_cps_ram,&rsnaps);CHKERRQ(ierr);
2227         revolve_create_offline(rfine,rsnaps);
2228         break;
2229       case REVOLVE_ONLINE:
2230         stack->stacksize = tjsch->max_cps_ram;
2231         ierr = PetscRevolveIntCast(tjsch->max_cps_ram,&rsnaps);CHKERRQ(ierr);
2232         revolve_create_online(rsnaps);
2233         break;
2234       case REVOLVE_MULTISTAGE:
2235         ierr = PetscRevolveIntCast(tjsch->total_steps,&rfine);CHKERRQ(ierr);
2236         ierr = PetscRevolveIntCast(tjsch->max_cps_ram,&rsnaps);CHKERRQ(ierr);
2237         ierr = PetscRevolveIntCast(tjsch->max_cps_ram+tjsch->max_cps_disk,&rsnaps2);CHKERRQ(ierr);
2238         revolve_create_multistage(rfine,rsnaps2,rsnaps);
2239         break;
2240       default:
2241         break;
2242     }
2243     ierr = PetscNew(&rctx);CHKERRQ(ierr);
2244     ierr = PetscRevolveIntCast(tjsch->max_cps_ram,&rsnaps);CHKERRQ(ierr);
2245     rctx->snaps_in       = rsnaps; /* for theta methods snaps_in=2*max_cps_ram */
2246     rctx->reverseonestep = PETSC_FALSE;
2247     rctx->check          = 0;
2248     rctx->oldcapo        = 0;
2249     rctx->capo           = 0;
2250     rctx->info           = 2;
2251     if (tjsch->stride > 1) {
2252       ierr = PetscRevolveIntCast(tjsch->stride,&rfine);CHKERRQ(ierr);
2253     } else {
2254       ierr = PetscRevolveIntCast(tjsch->total_steps,&rfine);CHKERRQ(ierr);
2255     }
2256     rctx->fine           = rfine;
2257     tjsch->rctx          = rctx;
2258     if (tjsch->stype == REVOLVE_ONLINE) rctx->fine = -1;
2259 #endif
2260   } else {
2261     if (tjsch->stype == TWO_LEVEL_NOREVOLVE) stack->stacksize = tjsch->stride-1; /* need tjsch->stride-1 at most */
2262     if (tjsch->stype == NONE) {
2263       if (fixedtimestep) stack->stacksize = stack->solution_only ? tjsch->total_steps : tjsch->total_steps-1;
2264       else { /* adaptive time step */
2265         /* if max_cps_ram is not specified, use maximal allowed number of steps for stack size */
2266         if (tjsch->max_cps_ram == -1) stack->stacksize = ts->max_steps < PETSC_MAX_INT ? ts->max_steps : 10000;
2267         tjsch->total_steps = stack->solution_only ? stack->stacksize : stack->stacksize+1; /* will be updated as time integration advances */
2268       }
2269     }
2270   }
2271 
2272   if ((tjsch->stype >= TWO_LEVEL_NOREVOLVE && tjsch->stype < REVOLVE_OFFLINE) || tjsch->stype == REVOLVE_MULTISTAGE) { /* these types need to use disk */
2273     ierr = TSTrajectorySetUp_Basic(tj,ts);CHKERRQ(ierr);
2274   }
2275 
2276   stack->stacksize = PetscMax(stack->stacksize,1);
2277   tjsch->recompute = PETSC_FALSE;
2278   ierr = StackInit(stack,stack->stacksize,numY);CHKERRQ(ierr);
2279   PetscFunctionReturn(0);
2280 }
2281 
2282 static PetscErrorCode TSTrajectoryReset_Memory(TSTrajectory tj)
2283 {
2284 #if defined (PETSC_HAVE_REVOLVE) || defined (PETSC_HAVE_CAMS)
2285   TJScheduler    *tjsch = (TJScheduler*)tj->data;
2286   PetscErrorCode ierr;
2287 #endif
2288 
2289   PetscFunctionBegin;
2290 #if defined(PETSC_HAVE_REVOLVE)
2291   if (tjsch->stype > TWO_LEVEL_NOREVOLVE) {
2292     revolve_reset();
2293     if (tjsch->stype == TWO_LEVEL_TWO_REVOLVE) {
2294       revolve2_reset();
2295       ierr = PetscFree(tjsch->diskstack.container);CHKERRQ(ierr);
2296     }
2297   }
2298   if (tjsch->stype > TWO_LEVEL_NOREVOLVE) {
2299     ierr = PetscFree(tjsch->rctx);CHKERRQ(ierr);
2300     ierr = PetscFree(tjsch->rctx2);CHKERRQ(ierr);
2301   }
2302 #endif
2303 #if defined(PETSC_HAVE_CAMS)
2304   if (tjsch->stype == CAMS_OFFLINE) {
2305     if (tjsch->stack.solution_only) offline_ca_destroy();
2306     else offline_ca_destroy();
2307     ierr = PetscFree(tjsch->actx);CHKERRQ(ierr);
2308   }
2309 #endif
2310   PetscFunctionReturn(0);
2311 }
2312 
2313 static PetscErrorCode TSTrajectoryDestroy_Memory(TSTrajectory tj)
2314 {
2315   TJScheduler    *tjsch = (TJScheduler*)tj->data;
2316   PetscErrorCode ierr;
2317 
2318   PetscFunctionBegin;
2319   ierr = StackDestroy(&tjsch->stack);CHKERRQ(ierr);
2320   ierr = PetscViewerDestroy(&tjsch->viewer);CHKERRQ(ierr);
2321   ierr = PetscObjectComposeFunction((PetscObject)tj,"TSTrajectorySetMaxCpsRAM_C",NULL);CHKERRQ(ierr);
2322   ierr = PetscObjectComposeFunction((PetscObject)tj,"TSTrajectorySetMaxCpsDisk_C",NULL);CHKERRQ(ierr);
2323   ierr = PetscObjectComposeFunction((PetscObject)tj,"TSTrajectorySetMaxUnitsRAM_C",NULL);CHKERRQ(ierr);
2324   ierr = PetscObjectComposeFunction((PetscObject)tj,"TSTrajectorySetMaxUnitsDisk_C",NULL);CHKERRQ(ierr);
2325   ierr = PetscFree(tjsch);CHKERRQ(ierr);
2326   PetscFunctionReturn(0);
2327 }
2328 
2329 /*MC
2330       TSTRAJECTORYMEMORY - Stores each solution of the ODE/ADE in memory
2331 
2332   Level: intermediate
2333 
2334 .seealso:  TSTrajectoryCreate(), TS, TSTrajectorySetType()
2335 
2336 M*/
2337 PETSC_EXTERN PetscErrorCode TSTrajectoryCreate_Memory(TSTrajectory tj,TS ts)
2338 {
2339   TJScheduler    *tjsch;
2340   PetscErrorCode ierr;
2341 
2342   PetscFunctionBegin;
2343   tj->ops->set            = TSTrajectorySet_Memory;
2344   tj->ops->get            = TSTrajectoryGet_Memory;
2345   tj->ops->setup          = TSTrajectorySetUp_Memory;
2346   tj->ops->setfromoptions = TSTrajectorySetFromOptions_Memory;
2347   tj->ops->reset          = TSTrajectoryReset_Memory;
2348   tj->ops->destroy        = TSTrajectoryDestroy_Memory;
2349 
2350   ierr = PetscNew(&tjsch);CHKERRQ(ierr);
2351   tjsch->stype        = NONE;
2352   tjsch->max_cps_ram  = -1; /* -1 indicates that it is not set */
2353   tjsch->max_cps_disk = -1; /* -1 indicates that it is not set */
2354   tjsch->stride       = 0; /* if not zero, two-level checkpointing will be used */
2355 #if defined(PETSC_HAVE_REVOLVE)
2356   tjsch->use_online   = PETSC_FALSE;
2357 #endif
2358   tjsch->save_stack   = PETSC_TRUE;
2359 
2360   tjsch->stack.solution_only = tj->solution_only;
2361   ierr = PetscViewerCreate(PetscObjectComm((PetscObject)tj),&tjsch->viewer);CHKERRQ(ierr);
2362   ierr = PetscViewerSetType(tjsch->viewer,PETSCVIEWERBINARY);CHKERRQ(ierr);
2363   ierr = PetscViewerPushFormat(tjsch->viewer,PETSC_VIEWER_NATIVE);CHKERRQ(ierr);
2364   ierr = PetscViewerFileSetMode(tjsch->viewer,FILE_MODE_WRITE);CHKERRQ(ierr);
2365 
2366   ierr = PetscObjectComposeFunction((PetscObject)tj,"TSTrajectorySetMaxCpsRAM_C",TSTrajectorySetMaxCpsRAM_Memory);CHKERRQ(ierr);
2367   ierr = PetscObjectComposeFunction((PetscObject)tj,"TSTrajectorySetMaxCpsDisk_C",TSTrajectorySetMaxCpsDisk_Memory);CHKERRQ(ierr);
2368   ierr = PetscObjectComposeFunction((PetscObject)tj,"TSTrajectorySetMaxUnitsRAM_C",TSTrajectorySetMaxUnitsRAM_Memory);CHKERRQ(ierr);
2369   ierr = PetscObjectComposeFunction((PetscObject)tj,"TSTrajectorySetMaxUnitsDisk_C",TSTrajectorySetMaxUnitsDisk_Memory);CHKERRQ(ierr);
2370   tj->data = tjsch;
2371   PetscFunctionReturn(0);
2372 }
2373