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