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