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