xref: /petsc/src/ts/trajectory/impls/memory/trajmemory.c (revision 00aea19f3afbf245e3a7fffa90018a8b125b1b19)
1 #include <petsc/private/tsimpl.h>        /*I "petscts.h"  I*/
2 #include <petscsys.h>
3 
4 #ifdef PETSC_HAVE_REVOLVE
5 #include <revolve_c.h>
6 #endif
7 
8 PetscLogEvent Disk_Write, Disk_Read;
9 
10 typedef enum {NONE,TWO_LEVEL_NOREVOLVE,TWO_LEVEL_REVOLVE,REVOLVE_OFFLINE,REVOLVE_ONLINE,REVOLVE_MULTISTAGE} SchedulerType;
11 
12 typedef struct _StackElement {
13   PetscInt  stepnum;
14   Vec       X;
15   Vec       *Y;
16   PetscReal time;
17   PetscReal timeprev;
18 } *StackElement;
19 
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 
32 typedef struct _Stack {
33   SchedulerType stype;
34   PetscBool     use_online;
35   PetscBool     recompute;
36   PetscBool     solution_only;
37   PetscBool     save_stack;
38   MPI_Comm      comm;
39   RevolveCTX    *rctx;
40   PetscInt      max_cps_ram;  /* maximum checkpoints in RAM */
41   PetscInt      max_cps_disk; /* maximum checkpoints on disk */
42   PetscInt      stride;
43   PetscInt      total_steps;  /* total number of steps */
44   PetscInt      numY;
45   PetscInt      stacksize;
46   PetscInt      top;          /* top of the stack */
47   StackElement  *stack;       /* container */
48 } Stack;
49 
50 #ifdef PETSC_HAVE_REVOLVE
51 static void printwhattodo(PetscInt whattodo,RevolveCTX *rctx,PetscInt shift)
52 {
53   switch(whattodo) {
54     case 1:
55       PetscPrintf(PETSC_COMM_WORLD,"\x1B[35mAdvance from %D to %D\033[0m\n",rctx->oldcapo+shift,rctx->capo+shift);
56       break;
57     case 2:
58       PetscPrintf(PETSC_COMM_WORLD,"\x1B[35mStore in checkpoint number %D (located in RAM)\033[0m\n",rctx->check);
59       break;
60     case 3:
61       PetscPrintf(PETSC_COMM_WORLD,"\x1B[35mFirst turn: Initialize adjoints and reverse first step\033[0m\n");
62       break;
63     case 4:
64       PetscPrintf(PETSC_COMM_WORLD,"\x1B[35mForward and reverse one step\033[0m\n");
65       break;
66     case 5:
67       PetscPrintf(PETSC_COMM_WORLD,"\x1B[35mRestore in checkpoint number %D (located in RAM)\033[0m\n",rctx->check);
68       break;
69     case 7:
70       PetscPrintf(PETSC_COMM_WORLD,"\x1B[35mStore in checkpoint number %D (located on disk)\033[0m\n",rctx->check);
71       break;
72     case 8:
73       PetscPrintf(PETSC_COMM_WORLD,"\x1B[35mRestore in checkpoint number %D (located on disk)\033[0m\n",rctx->check);
74       break;
75     case -1:
76       PetscPrintf(PETSC_COMM_WORLD,"\x1B[35mError!");
77       break;
78   }
79 }
80 #endif
81 
82 #undef __FUNCT__
83 #define __FUNCT__ "StackCreate"
84 static PetscErrorCode StackCreate(MPI_Comm comm,Stack *s,PetscInt size,PetscInt ny)
85 {
86   PetscErrorCode ierr;
87 
88   PetscFunctionBegin;
89   s->top  = -1;
90   s->comm = comm;
91   s->numY = ny;
92 
93   ierr = PetscMalloc1(size*sizeof(StackElement),&s->stack);CHKERRQ(ierr);
94   PetscFunctionReturn(0);
95 }
96 
97 #undef __FUNCT__
98 #define __FUNCT__ "StackDestroy"
99 static PetscErrorCode StackDestroy(Stack *s)
100 {
101   PetscInt       i;
102   PetscErrorCode ierr;
103 
104   PetscFunctionBegin;
105   if (s->top>-1) {
106     for (i=0;i<=s->top;i++) {
107       ierr = VecDestroy(&s->stack[i]->X);CHKERRQ(ierr);
108       if (!s->solution_only) {
109         ierr = VecDestroyVecs(s->numY,&s->stack[i]->Y);CHKERRQ(ierr);
110       }
111       ierr = PetscFree(s->stack[i]);CHKERRQ(ierr);
112     }
113   }
114   ierr = PetscFree(s->stack);CHKERRQ(ierr);
115   if (s->stype) {
116     ierr = PetscFree(s->rctx);CHKERRQ(ierr);
117   }
118   ierr = PetscFree(s);CHKERRQ(ierr);
119   PetscFunctionReturn(0);
120 }
121 
122 #undef __FUNCT__
123 #define __FUNCT__ "StackPush"
124 static PetscErrorCode StackPush(Stack *s,StackElement e)
125 {
126   PetscFunctionBegin;
127   if (s->top+1 >= s->stacksize) SETERRQ1(s->comm,PETSC_ERR_MEMC,"Maximum stack size (%D) exceeded",s->stacksize);
128   s->stack[++s->top] = e;
129   PetscFunctionReturn(0);
130 }
131 
132 #undef __FUNCT__
133 #define __FUNCT__ "StackPop"
134 static PetscErrorCode StackPop(Stack *s,StackElement *e)
135 {
136   PetscFunctionBegin;
137   if (s->top == -1) SETERRQ(s->comm,PETSC_ERR_MEMC,"Empty stack");
138   *e = s->stack[s->top--];
139   PetscFunctionReturn(0);
140 }
141 
142 #undef __FUNCT__
143 #define __FUNCT__ "StackTop"
144 static PetscErrorCode StackTop(Stack *s,StackElement *e)
145 {
146   PetscFunctionBegin;
147   *e = s->stack[s->top];
148   PetscFunctionReturn(0);
149 }
150 
151 #undef __FUNCT__
152 #define __FUNCT__ "StackFind"
153 static PetscErrorCode StackFind(Stack *s,StackElement *e,PetscInt index)
154 {
155   PetscFunctionBegin;
156   *e = s->stack[index];
157   PetscFunctionReturn(0);
158 }
159 
160 #undef __FUNCT__
161 #define __FUNCT__ "OutputBIN"
162 static PetscErrorCode OutputBIN(const char *filename,PetscViewer *viewer)
163 {
164   PetscErrorCode ierr;
165 
166   PetscFunctionBegin;
167   ierr = PetscViewerCreate(PETSC_COMM_WORLD,viewer);CHKERRQ(ierr);
168   ierr = PetscViewerSetType(*viewer,PETSCVIEWERBINARY);CHKERRQ(ierr);
169   ierr = PetscViewerFileSetMode(*viewer,FILE_MODE_WRITE);CHKERRQ(ierr);
170   ierr = PetscViewerFileSetName(*viewer,filename);CHKERRQ(ierr);
171   PetscFunctionReturn(0);
172 }
173 
174 #undef __FUNCT__
175 #define __FUNCT__ "WriteToDisk"
176 static PetscErrorCode WriteToDisk(PetscInt stepnum,PetscReal time,PetscReal timeprev,Vec X,Vec *Y,PetscInt numY,PetscBool solution_only,PetscViewer viewer)
177 {
178   PetscInt       i;
179   PetscErrorCode ierr;
180 
181   PetscFunctionBegin;
182   ierr = PetscViewerBinaryWrite(viewer,&stepnum,1,PETSC_INT,PETSC_FALSE);CHKERRQ(ierr);
183   ierr = VecView(X,viewer);CHKERRQ(ierr);
184   for (i=0;!solution_only && i<numY;i++) {
185     ierr = VecView(Y[i],viewer);CHKERRQ(ierr);
186   }
187   ierr = PetscViewerBinaryWrite(viewer,&time,1,PETSC_REAL,PETSC_FALSE);CHKERRQ(ierr);
188   ierr = PetscViewerBinaryWrite(viewer,&timeprev,1,PETSC_REAL,PETSC_FALSE);CHKERRQ(ierr);
189   PetscFunctionReturn(0);
190 }
191 
192 #undef __FUNCT__
193 #define __FUNCT__ "ReadFromDisk"
194 static PetscErrorCode ReadFromDisk(PetscInt *stepnum,PetscReal *time,PetscReal *timeprev,Vec X,Vec *Y,PetscInt numY,PetscBool solution_only,PetscViewer viewer)
195 {
196   PetscInt       i;
197   PetscErrorCode ierr;
198 
199   PetscFunctionBegin;
200   ierr = PetscViewerBinaryRead(viewer,stepnum,1,NULL,PETSC_INT);CHKERRQ(ierr);
201   ierr = VecLoad(X,viewer);CHKERRQ(ierr);
202   for (i=0;!solution_only && i<numY;i++) {
203     ierr = VecLoad(Y[i],viewer);CHKERRQ(ierr);
204   }
205   ierr = PetscViewerBinaryRead(viewer,time,1,NULL,PETSC_REAL);CHKERRQ(ierr);
206   ierr = PetscViewerBinaryRead(viewer,timeprev,1,NULL,PETSC_REAL);CHKERRQ(ierr);
207   PetscFunctionReturn(0);
208 }
209 
210 #undef __FUNCT__
211 #define __FUNCT__ "StackDumpAll"
212 static PetscErrorCode StackDumpAll(TS ts,Stack *s,PetscInt id)
213 {
214   Vec            *Y;
215   PetscInt       i;
216   StackElement   e;
217   PetscViewer    viewer;
218   char           filename[PETSC_MAX_PATH_LEN];
219   PetscErrorCode ierr;
220 
221   PetscFunctionBegin;
222   if (id == 1) {
223     PetscMPIInt rank;
224     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)ts),&rank);CHKERRQ(ierr);
225     if (!rank) {
226       ierr = PetscRMTree("SA-data");CHKERRQ(ierr);
227       ierr = PetscMkdir("SA-data");CHKERRQ(ierr);
228     }
229   }
230   ierr = PetscSNPrintf(filename,sizeof(filename),"SA-data/SA-STACK%06d.bin",id);CHKERRQ(ierr);
231   ierr = OutputBIN(filename,&viewer);CHKERRQ(ierr);
232   for (i=0;i<s->stacksize;i++) {
233     e = s->stack[i];
234     ierr = WriteToDisk(e->stepnum,e->time,e->timeprev,e->X,e->Y,s->numY,s->solution_only,viewer);CHKERRQ(ierr);
235   }
236   /* 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 */
237   ierr = TSGetStages(ts,&s->numY,&Y);CHKERRQ(ierr);
238   ierr = WriteToDisk(ts->total_steps,ts->ptime,ts->ptime_prev,ts->vec_sol,Y,s->numY,s->solution_only,viewer);CHKERRQ(ierr);
239   for (i=0;i<s->stacksize;i++) {
240     ierr = StackPop(s,&e);CHKERRQ(ierr);
241     ierr = VecDestroy(&e->X);CHKERRQ(ierr);
242     if (!s->solution_only) {
243       ierr = VecDestroyVecs(s->numY,&e->Y);CHKERRQ(ierr);
244     }
245     ierr = PetscFree(e);CHKERRQ(ierr);
246   }
247   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
248   PetscFunctionReturn(0);
249 }
250 
251 #undef __FUNCT__
252 #define __FUNCT__ "StackLoadAll"
253 static PetscErrorCode StackLoadAll(TS ts,Stack *s,PetscInt id)
254 {
255   Vec            *Y;
256   PetscInt       i;
257   StackElement   e;
258   PetscViewer    viewer;
259   char           filename[PETSC_MAX_PATH_LEN];
260   PetscErrorCode ierr;
261 
262   PetscFunctionBegin;
263   ierr = PetscSNPrintf(filename,sizeof filename,"SA-data/SA-STACK%06d.bin",id);CHKERRQ(ierr);
264   ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);CHKERRQ(ierr);
265   for (i=0;i<s->stacksize;i++) {
266     ierr = PetscCalloc1(1,&e);
267     ierr = TSGetStages(ts,&s->numY,&Y);CHKERRQ(ierr);
268     ierr = VecDuplicate(Y[0],&e->X);CHKERRQ(ierr);
269     if (!s->solution_only && s->numY>0) {
270       ierr = VecDuplicateVecs(Y[0],s->numY,&e->Y);CHKERRQ(ierr);
271     }
272     ierr = StackPush(s,e);CHKERRQ(ierr);
273   }
274   for (i=0;i<s->stacksize;i++) {
275     e = s->stack[i];
276     ierr = ReadFromDisk(&e->stepnum,&e->time,&e->timeprev,e->X,e->Y,s->numY,s->solution_only,viewer);CHKERRQ(ierr);
277   }
278   /* load the last step into TS */
279   ierr = TSGetStages(ts,&s->numY,&Y);CHKERRQ(ierr);
280   ierr = ReadFromDisk(&ts->total_steps,&ts->ptime,&ts->ptime_prev,ts->vec_sol,Y,s->numY,s->solution_only,viewer);CHKERRQ(ierr);
281   ierr = TSSetTimeStep(ts,ts->ptime-ts->ptime_prev);CHKERRQ(ierr);
282   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
283   PetscFunctionReturn(0);
284 }
285 
286 #undef __FUNCT__
287 #define __FUNCT__ "DumpSingle"
288 static PetscErrorCode DumpSingle(TS ts,Stack *s,PetscInt id)
289 {
290   Vec            *Y;
291   PetscInt       stepnum;
292   PetscViewer    viewer;
293   char           filename[PETSC_MAX_PATH_LEN];
294   PetscErrorCode ierr;
295 
296   PetscFunctionBegin;
297   ierr = TSGetTotalSteps(ts,&stepnum);CHKERRQ(ierr);
298   if (id == 0) {
299     PetscMPIInt rank;
300     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)ts),&rank);CHKERRQ(ierr);
301     if (!rank) {
302       ierr = PetscRMTree("SA-data");CHKERRQ(ierr);
303       ierr = PetscMkdir("SA-data");CHKERRQ(ierr);
304     }
305   }
306   ierr = PetscSNPrintf(filename,sizeof(filename),"SA-data/SA-CPS%06d.bin",id);CHKERRQ(ierr);
307   ierr = OutputBIN(filename,&viewer);CHKERRQ(ierr);
308 
309   ierr = PetscLogEventBegin(Disk_Write,ts,0,0,0);CHKERRQ(ierr);
310   ierr = TSGetStages(ts,&s->numY,&Y);CHKERRQ(ierr);
311   ierr = WriteToDisk(stepnum,ts->ptime,ts->ptime_prev,ts->vec_sol,Y,s->numY,s->solution_only,viewer);CHKERRQ(ierr);
312   ierr = PetscLogEventEnd(Disk_Write,ts,0,0,0);CHKERRQ(ierr);
313 
314   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
315   PetscFunctionReturn(0);
316 }
317 
318 #undef __FUNCT__
319 #define __FUNCT__ "LoadSingle"
320 static PetscErrorCode LoadSingle(TS ts,Stack *s,PetscInt id)
321 {
322   Vec            *Y;
323   PetscViewer    viewer;
324   char           filename[PETSC_MAX_PATH_LEN];
325   PetscErrorCode ierr;
326 
327   PetscFunctionBegin;
328   ierr = PetscSNPrintf(filename,sizeof filename,"SA-data/SA-CPS%06d.bin",id);CHKERRQ(ierr);
329   ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);CHKERRQ(ierr);
330 
331   ierr = PetscLogEventBegin(Disk_Read,ts,0,0,0);CHKERRQ(ierr);
332   ierr = TSGetStages(ts,&s->numY,&Y);CHKERRQ(ierr);
333   ierr = ReadFromDisk(&ts->total_steps,&ts->ptime,&ts->ptime_prev,ts->vec_sol,Y,s->numY,s->solution_only,viewer);CHKERRQ(ierr);
334   ierr = PetscLogEventEnd(Disk_Read,ts,0,0,0);CHKERRQ(ierr);
335 
336   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
337   PetscFunctionReturn(0);
338 }
339 
340 #undef __FUNCT__
341 #define __FUNCT__ "TSTrajectorySetStride_Memory"
342 PetscErrorCode TSTrajectorySetStride_Memory(TSTrajectory tj,TS ts,PetscInt stride)
343 {
344   Stack *s = (Stack*)tj->data;
345 
346   PetscFunctionBegin;
347   s->stride = stride;
348   PetscFunctionReturn(0);
349 }
350 
351 #undef __FUNCT__
352 #define __FUNCT__ "TSTrajectorySetMaxCpsRAM_Memory"
353 PetscErrorCode TSTrajectorySetMaxCpsRAM_Memory(TSTrajectory tj,TS ts,PetscInt max_cps_ram)
354 {
355   Stack *s = (Stack*)tj->data;
356 
357   PetscFunctionBegin;
358   s->max_cps_ram = max_cps_ram;
359   PetscFunctionReturn(0);
360 }
361 
362 #undef __FUNCT__
363 #define __FUNCT__ "TSTrajectorySetMaxCpsDisk_Memory"
364 PetscErrorCode TSTrajectorySetMaxCpsDisk_Memory(TSTrajectory tj,TS ts,PetscInt max_cps_disk)
365 {
366   Stack *s = (Stack*)tj->data;
367 
368   PetscFunctionBegin;
369   s->max_cps_disk = max_cps_disk;
370   PetscFunctionReturn(0);
371 }
372 
373 #undef __FUNCT__
374 #define __FUNCT__ "TSTrajectorySetRevolveOnline"
375 PetscErrorCode TSTrajectorySetRevolveOnline(TSTrajectory tj,PetscBool use_online)
376 {
377   Stack *s = (Stack*)tj->data;
378 
379   PetscFunctionBegin;
380   s->use_online = use_online;
381   PetscFunctionReturn(0);
382 }
383 
384 #undef __FUNCT__
385 #define __FUNCT__ "TSTrajectorySetSaveStack"
386 PetscErrorCode TSTrajectorySetSaveStack(TSTrajectory tj,PetscBool save_stack)
387 {
388   Stack *s = (Stack*)tj->data;
389 
390   PetscFunctionBegin;
391   s->save_stack = save_stack;
392   PetscFunctionReturn(0);
393 }
394 
395 #undef __FUNCT__
396 #define __FUNCT__ "TSTrajectorySetSolutionOnly"
397 PetscErrorCode TSTrajectorySetSolutionOnly(TSTrajectory tj,PetscBool solution_only)
398 {
399   Stack *s = (Stack*)tj->data;
400 
401   PetscFunctionBegin;
402   s->solution_only = solution_only;
403   PetscFunctionReturn(0);
404 }
405 
406 #undef __FUNCT__
407 #define __FUNCT__ "TSTrajectorySetFromOptions_Memory"
408 PetscErrorCode TSTrajectorySetFromOptions_Memory(PetscOptions *PetscOptionsObject,TSTrajectory tj)
409 {
410   Stack          *s = (Stack*)tj->data;
411   PetscErrorCode ierr;
412 
413   PetscFunctionBegin;
414   ierr = PetscOptionsHead(PetscOptionsObject,"Memory based TS trajectory options");CHKERRQ(ierr);
415   {
416     ierr = PetscOptionsInt("-tstrajectory_max_cps_ram","Maximum number of checkpoints in RAM","TSTrajectorySetMaxCpsRAM_Memory",s->max_cps_ram,&s->max_cps_ram,NULL);CHKERRQ(ierr);
417     ierr = PetscOptionsInt("-tstrajectory_max_cps_disk","Maximum number of checkpoints on disk","TSTrajectorySetMaxCpsDisk_Memory",s->max_cps_disk,&s->max_cps_disk,NULL);CHKERRQ(ierr);
418     ierr = PetscOptionsInt("-tstrajectory_stride","Stride to save checkpoints to file","TSTrajectorySetStride_Memory",s->stride,&s->stride,NULL);CHKERRQ(ierr);
419     ierr = PetscOptionsBool("-tstrajectory_revolve_online","Trick TS trajectory into using online mode of revolve","TSTrajectorySetRevolveOnline",s->use_online,&s->use_online,NULL);CHKERRQ(ierr);
420     ierr = PetscOptionsBool("-tstrajectory_save_stack","Save all stack to disk","TSTrajectorySetSaveStack",s->save_stack,&s->save_stack,NULL);CHKERRQ(ierr);
421     ierr = PetscOptionsBool("-tstrajectory_solution_only","Checkpoint solution only","TSTrajectorySetSolutionOnly",s->solution_only,&s->solution_only,NULL);CHKERRQ(ierr);
422   }
423   ierr = PetscOptionsTail();CHKERRQ(ierr);
424   PetscFunctionReturn(0);
425 }
426 
427 #undef __FUNCT__
428 #define __FUNCT__ "TSTrajectorySetUp_Memory"
429 PetscErrorCode TSTrajectorySetUp_Memory(TSTrajectory tj,TS ts)
430 {
431   Stack          *s = (Stack*)tj->data;
432   RevolveCTX     *rctx;
433   PetscInt       numY;
434   PetscBool      flg;
435   PetscErrorCode ierr;
436 
437   PetscFunctionBegin;
438   PetscStrcmp(((PetscObject)ts->adapt)->type_name,TSADAPTNONE,&flg);
439   if (flg) { /* fixed time step */
440     s->total_steps = PetscMin(ts->max_steps,(PetscInt)(ceil(ts->max_time/ts->time_step)));
441   }
442   if (s->max_cps_ram > 1) s->stacksize = s->max_cps_ram;
443   if (s->stride > 1) { /* two level mode works for both fixed time step and adaptive time step */
444     if (s->max_cps_ram > 1 && s->max_cps_ram < s->stride-1) { /* use revolve_offline for each stride */
445       s->stype = TWO_LEVEL_REVOLVE;
446     }else { /* checkpoint all for each stride */
447       s->stype     = TWO_LEVEL_NOREVOLVE;
448       s->stacksize = s->solution_only ? s->stride : s->stride-1;
449     }
450   } else {
451     if (flg) { /* fixed time step */
452       if (s->max_cps_ram >= s->total_steps-1 || s->max_cps_ram < 1) { /* checkpoint all */
453         s->stype     = NONE;
454         s->stacksize = s->solution_only ? s->total_steps : s->total_steps-1;
455       } else {
456         if (s->max_cps_disk > 1) { /* disk can be used */
457           s->stype = REVOLVE_MULTISTAGE;
458         } else { /* memory only */
459           s->stype = REVOLVE_OFFLINE;
460         }
461       }
462     } else { /* adaptive time step */
463       s->stype = REVOLVE_ONLINE;
464     }
465     if (s->use_online) { /* trick into online */
466       s->stype     = REVOLVE_ONLINE;
467       s->stacksize = s->max_cps_ram;
468     }
469   }
470 
471   if (s->stype > TWO_LEVEL_NOREVOLVE) {
472 #ifndef PETSC_HAVE_REVOLVE
473     SETERRQ(s->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.");
474 #else
475     if (s->stype == TWO_LEVEL_REVOLVE) revolve_create_offline(s->stride,s->max_cps_ram);
476     else if (s->stype == REVOLVE_OFFLINE) revolve_create_offline(s->total_steps,s->max_cps_ram);
477     else if (s->stype == REVOLVE_ONLINE) revolve_create_online(s->max_cps_ram);
478     else if (s->stype ==REVOLVE_MULTISTAGE) revolve_create_multistage(s->total_steps,s->max_cps_ram+s->max_cps_disk,s->max_cps_ram);
479 
480     ierr = PetscCalloc1(1,&rctx);CHKERRQ(ierr);
481     rctx->snaps_in       = s->max_cps_ram; /* for theta methods snaps_in=2*max_cps_ram */
482     rctx->reverseonestep = PETSC_FALSE;
483     rctx->check          = 0;
484     rctx->oldcapo        = 0;
485     rctx->capo           = 0;
486     rctx->info           = 2;
487     rctx->fine           = (s->stride > 1) ? s->stride : s->total_steps;
488 
489     s->rctx      = rctx;
490     if (s->stype == REVOLVE_ONLINE) rctx->fine = -1;
491 #endif
492   }
493 
494   s->recompute = PETSC_FALSE;
495 
496   ierr = TSGetStages(ts,&numY,PETSC_IGNORE);CHKERRQ(ierr);
497   ierr = StackCreate(PetscObjectComm((PetscObject)ts),s,s->stacksize,numY);CHKERRQ(ierr);
498   PetscFunctionReturn(0);
499 }
500 
501 #undef __FUNCT__
502 #define __FUNCT__ "ApplyRevolve"
503 static PetscErrorCode ApplyRevolve(TS ts,Stack *s,PetscInt stepnum,PetscInt localstepnum,PetscInt *whattodo)
504 {
505 #ifdef PETSC_HAVE_REVOLVE
506   PetscInt       shift;
507 #endif
508   RevolveCTX     *rctx;
509   PetscErrorCode ierr;
510 
511   PetscFunctionBegin;
512 #ifdef PETSC_HAVE_REVOLVE
513     rctx = s->rctx;
514     if (rctx->reverseonestep && stepnum==s->total_steps) { /* intermediate information is ready inside TS, this happens at last time step */
515       //rctx->reverseonestep = PETSC_FALSE;
516       PetscFunctionReturn(0);
517     }
518     if (rctx->stepsleft != 0) { /* advance the solution without checkpointing anything as Revolve requires */
519       rctx->stepsleft--;
520       PetscFunctionReturn(0);
521     }
522 
523     /* let Revolve determine what to do next */
524     shift         = stepnum-localstepnum;
525     rctx->capo    = localstepnum;
526     rctx->oldcapo = rctx->capo;
527     *whattodo = revolve_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where);
528     if (s->stype == REVOLVE_ONLINE && *whattodo ==7) *whattodo = 2;
529     printwhattodo(*whattodo,rctx,shift);
530     if (*whattodo == -1) SETERRQ(s->comm,PETSC_ERR_LIB,"Error in the Revolve library");
531     if (*whattodo == 1) { /* advance some time steps */
532       if (s->stype == REVOLVE_ONLINE && rctx->capo >= s->total_steps-1) {
533         revolve_turn(s->total_steps,&rctx->capo,&rctx->fine);
534         printwhattodo(*whattodo,rctx,shift);
535       }
536       rctx->stepsleft = rctx->capo-rctx->oldcapo-1;
537     }
538     if (*whattodo == 3 || *whattodo == 4) { /* ready for a reverse step */
539       rctx->reverseonestep = PETSC_TRUE;
540     }
541     if (*whattodo == 5) { /* restore a checkpoint and ask Revolve what to do next */
542       rctx->oldcapo = rctx->capo;
543       *whattodo = revolve_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where); /* must return 1*/
544       printwhattodo(*whattodo,rctx,shift);
545       rctx->stepsleft = rctx->capo-rctx->oldcapo;
546     }
547     if (*whattodo == 7) { /* save the checkpoint to disk */
548       ierr = DumpSingle(ts,s,rctx->check);CHKERRQ(ierr);
549       rctx->oldcapo = rctx->capo;
550       *whattodo = revolve_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where); /* must return 1*/
551       printwhattodo(*whattodo,rctx,shift);
552       rctx->stepsleft = rctx->capo-rctx->oldcapo-1;
553     }
554     if (*whattodo != 2) {
555       PetscFunctionReturn(0);
556     } else { /* store a checkpoint to RAM and ask Revolve how many time steps to advance next */
557       rctx->oldcapo = rctx->capo;
558       *whattodo = revolve_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where); /* must return 1*/
559       printwhattodo(*whattodo,rctx,shift);
560       if (s->stype == REVOLVE_ONLINE && rctx->capo >= s->total_steps-1) {
561         revolve_turn(s->total_steps,&rctx->capo,&rctx->fine);
562         printwhattodo(*whattodo,rctx,shift);
563       }
564       rctx->stepsleft = rctx->capo-rctx->oldcapo-1;
565     }
566 #endif
567   PetscFunctionReturn(0);
568 }
569 
570 #undef __FUNCT__
571 #define __FUNCT__ "ElementCreate"
572 static PetscErrorCode ElementCreate(TS ts,Stack *s,StackElement *e,PetscInt stepnum,PetscReal time,Vec X)
573 {
574   Vec            *Y;
575   PetscInt       i;
576   PetscReal      timeprev;
577   PetscErrorCode ierr;
578 
579   PetscFunctionBegin;
580   ierr = PetscCalloc1(1,e);CHKERRQ(ierr);
581   ierr = VecDuplicate(X,&(*e)->X);CHKERRQ(ierr);
582   ierr = VecCopy(X,(*e)->X);CHKERRQ(ierr);
583   if (s->numY > 0 && !s->solution_only) {
584     ierr = TSGetStages(ts,&s->numY,&Y);CHKERRQ(ierr);
585     ierr = VecDuplicateVecs(Y[0],s->numY,&(*e)->Y);CHKERRQ(ierr);
586     for (i=0;i<s->numY;i++) {
587       ierr = VecCopy(Y[i],(*e)->Y[i]);CHKERRQ(ierr);
588     }
589   }
590   (*e)->stepnum = stepnum;
591   (*e)->time    = time;
592   /* for consistency */
593   if (stepnum == 0) {
594     (*e)->timeprev = (*e)->time - ts->time_step;
595   } else {
596     ierr = TSGetPrevTime(ts,&timeprev);CHKERRQ(ierr);
597     (*e)->timeprev = timeprev;
598   }
599   PetscFunctionReturn(0);
600 }
601 
602 #undef __FUNCT__
603 #define __FUNCT__ "ElementDestroy"
604 static PetscErrorCode ElementDestroy(Stack *s,StackElement e)
605 {
606   PetscErrorCode ierr;
607 
608   PetscFunctionBegin;
609   ierr = VecDestroy(&e->X);CHKERRQ(ierr);
610   if (!s->solution_only) {
611     ierr = VecDestroyVecs(s->numY,&e->Y);CHKERRQ(ierr);
612   }
613   ierr = PetscFree(e);CHKERRQ(ierr);
614   PetscFunctionReturn(0);
615 }
616 
617 #undef __FUNCT__
618 #define __FUNCT__ "UpdateTS"
619 static PetscErrorCode UpdateTS(TS ts,Stack *s,StackElement e)
620 {
621   Vec            *Y;
622   PetscInt       i;
623   PetscErrorCode ierr;
624 
625   PetscFunctionBegin;
626   ierr = VecCopy(e->X,ts->vec_sol);CHKERRQ(ierr);
627   if (!s->solution_only) {
628     ierr = TSGetStages(ts,&s->numY,&Y);CHKERRQ(ierr);
629     for (i=0;i<s->numY;i++) {
630       ierr = VecCopy(e->Y[i],Y[i]);CHKERRQ(ierr);
631     }
632   }
633   ierr = TSSetTimeStep(ts,e->timeprev-e->time);CHKERRQ(ierr); /* stepsize will be negative */
634   PetscFunctionReturn(0);
635 }
636 
637 #undef __FUNCT__
638 #define __FUNCT__ "ReCompute"
639 static PetscErrorCode ReCompute(TS ts,Stack *s,StackElement e,PetscInt stepnum)
640 {
641   PetscInt       i;
642   PetscErrorCode ierr;
643 
644   PetscFunctionBegin;
645   /* reset ts context */
646   ts->steps      = e->stepnum; /* global stepnum */
647   ts->ptime      = e->time;
648   ts->ptime_prev = e->timeprev;
649   for (i=ts->steps;i<stepnum;i++) { /* assume fixed step size */
650     ierr = TSTrajectorySet(ts->trajectory,ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
651     ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
652     ierr = TSStep(ts);CHKERRQ(ierr);
653     if (ts->event) {
654       ierr = TSEventMonitor(ts);CHKERRQ(ierr);
655     }
656     if (!ts->steprollback) {
657       ierr = TSPostStep(ts);CHKERRQ(ierr);
658     }
659   }
660   PetscFunctionReturn(0);
661 }
662 
663 #undef __FUNCT__
664 #define __FUNCT__ "SetTrajN"
665 static PetscErrorCode SetTrajN(TS ts,Stack *s,PetscInt stepnum,PetscReal time,Vec X)
666 {
667   StackElement   e;
668   PetscErrorCode ierr;
669 
670   PetscFunctionBegin;
671   if (s->solution_only) {
672     /* skip the last two steps of each stride or the whole interval */
673     if (stepnum >= s->total_steps-1 || s->recompute) PetscFunctionReturn(0); //?
674   } else {
675     /* skip the first and the last steps of each stride or the whole interval */
676     if (stepnum == 0 || stepnum == s->total_steps) PetscFunctionReturn(0);
677   }
678   if (stepnum < s->top) SETERRQ(s->comm,PETSC_ERR_MEMC,"Illegal modification of a non-top stack element");
679   ierr = ElementCreate(ts,s,&e,stepnum,time,X);CHKERRQ(ierr);
680   ierr = StackPush(s,e);CHKERRQ(ierr);
681   PetscFunctionReturn(0);
682 }
683 
684 #undef __FUNCT__
685 #define __FUNCT__ "GetTrajN"
686 static PetscErrorCode GetTrajN(TS ts,Stack *s,PetscInt stepnum)
687 {
688   PetscInt       steps;
689   PetscReal      stepsize;
690   StackElement   e;
691   PetscErrorCode ierr;
692 
693   PetscFunctionBegin;
694   if (stepnum == s->total_steps) {
695     ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr);
696     ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr);
697     PetscFunctionReturn(0);
698   }
699   /* restore a checkpoint */
700   ierr = StackTop(s,&e);CHKERRQ(ierr);
701   ierr = UpdateTS(ts,s,e);CHKERRQ(ierr);
702   if (s->solution_only) {/* recompute one step */
703     steps = ts->steps;
704     s->recompute = PETSC_TRUE;
705     ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr);
706     ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr);
707     ierr = ReCompute(ts,s,e,stepnum);CHKERRQ(ierr);
708     ts->steps = steps;
709     //ts->total_steps = stepnum;
710     ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr);
711     ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr);
712   }
713   ierr = StackPop(s,&e);CHKERRQ(ierr);
714   ierr = ElementDestroy(s,e);CHKERRQ(ierr);
715   PetscFunctionReturn(0);
716 }
717 
718 #undef __FUNCT__
719 #define __FUNCT__ "SetTrajROF"
720 static PetscErrorCode SetTrajROF(TS ts,Stack *s,PetscInt stepnum,PetscReal time,Vec X)
721 {
722   PetscInt       whattodo;
723   StackElement   e;
724   PetscErrorCode ierr;
725 
726   PetscFunctionBegin;
727   ierr = ApplyRevolve(ts,s,stepnum,stepnum,&whattodo);CHKERRQ(ierr);
728   if (whattodo == 2){
729     if (stepnum < s->top) SETERRQ(s->comm,PETSC_ERR_MEMC,"Illegal modification of a non-top stack element");
730     ierr = ElementCreate(ts,s,&e,stepnum,time,X);CHKERRQ(ierr);
731     ierr = StackPush(s,e);CHKERRQ(ierr);
732   }
733   PetscFunctionReturn(0);
734 }
735 
736 #undef __FUNCT__
737 #define __FUNCT__ "GetTrajROF"
738 static PetscErrorCode GetTrajROF(TS ts,Stack *s,PetscInt stepnum)
739 {
740 #ifdef PETSC_HAVE_REVOLVE
741   PetscInt       whattodo,shift;
742 #endif
743   PetscInt       steps;
744   PetscReal      stepsize;
745   StackElement   e;
746   PetscErrorCode ierr;
747 
748   PetscFunctionBegin;
749 #ifdef PETSC_HAVE_REVOLVE
750   if ( s->rctx->reverseonestep) { /* ready for the reverse step */
751     ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr);
752     ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr);
753     s->rctx->reverseonestep = PETSC_FALSE;
754     PetscFunctionReturn(0);
755   }
756 #endif
757   if ((!s->solution_only && stepnum == 0) || stepnum == s->total_steps) {
758     ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr);
759     ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr);
760     PetscFunctionReturn(0);
761   }
762 
763   if (s->solution_only) {
764 #ifdef PETSC_HAVE_REVOLVE
765     s->rctx->capo = stepnum;
766     shift = 0;
767     whattodo = revolve_action(&s->rctx->check,&s->rctx->capo,&s->rctx->fine,s->rctx->snaps_in,&s->rctx->info,&s->rctx->where);
768     printwhattodo(whattodo,s->rctx,shift);
769 #endif
770   }
771   /* restore a checkpoint */
772   ierr = StackTop(s,&e);CHKERRQ(ierr);
773   if (e && e->stepnum >= stepnum) {
774     SETERRQ2(s->comm,PETSC_ERR_ARG_OUTOFRANGE,"The current step no. is %D, but the step number at top of the stack is %D",stepnum,e->stepnum);
775   }
776   ierr = UpdateTS(ts,s,e);CHKERRQ(ierr);
777   if (s->solution_only || (!s->solution_only && e->stepnum<stepnum)) { /* must recompute */
778     steps = ts->steps;
779     s->recompute = PETSC_TRUE;
780     ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr);
781     ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr);
782     ierr = ReCompute(ts,s,e,stepnum);CHKERRQ(ierr);
783     ts->steps = steps;
784     //ts->total_steps = stepnum;
785     ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr);
786     ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr);
787     #ifdef PETSC_HAVE_REVOLVE
788       s->rctx->reverseonestep = PETSC_FALSE;
789     #endif
790   } else {
791     ierr = ElementDestroy(s,e);CHKERRQ(ierr);
792   }
793   PetscFunctionReturn(0);
794 }
795 
796 #undef __FUNCT__
797 #define __FUNCT__ "SetTrajRON"
798 static PetscErrorCode SetTrajRON(TS ts,Stack *s,PetscInt stepnum,PetscReal time,Vec X)
799 {
800   PetscInt       whattodo;
801   PetscReal      timeprev;
802   StackElement   e;
803   RevolveCTX     *rctx = s->rctx;
804   PetscErrorCode ierr;
805 
806   PetscFunctionBegin;
807   ierr = ApplyRevolve(ts,s,stepnum,stepnum,&whattodo);CHKERRQ(ierr);
808   if (whattodo == 2){
809     if (rctx->check != s->top+1) { /* overwrite some non-top checkpoint in the stack*/
810       ierr = StackFind(s,&e,rctx->check);CHKERRQ(ierr);
811       ierr = VecCopy(X,e->X);CHKERRQ(ierr);
812       e->stepnum  = stepnum;
813       e->time     = time;
814       ierr        = TSGetPrevTime(ts,&timeprev);CHKERRQ(ierr);
815       e->timeprev = timeprev;
816     } else {
817       if (stepnum < s->top) SETERRQ(s->comm,PETSC_ERR_MEMC,"Illegal modification of a non-top stack element");
818       ierr = ElementCreate(ts,s,&e,stepnum,time,X);CHKERRQ(ierr);
819       ierr = StackPush(s,e);CHKERRQ(ierr);
820     }
821   }
822   PetscFunctionReturn(0);
823 }
824 
825 #undef __FUNCT__
826 #define __FUNCT__ "GetTrajRON"
827 static PetscErrorCode GetTrajRON(TS ts,Stack *s,PetscInt stepnum)
828 {
829 #ifdef PETSC_HAVE_REVOLVE
830   PetscInt       whattodo,shift;
831 #endif
832   PetscInt       steps;
833   PetscReal      stepsize;
834   StackElement   e;
835   PetscErrorCode ierr;
836 
837   PetscFunctionBegin;
838 #ifdef PETSC_HAVE_REVOLVE
839   if ( s->rctx->reverseonestep) { /* ready for the reverse step */
840     ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr);
841     ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr);
842     s->rctx->reverseonestep = PETSC_FALSE;
843     PetscFunctionReturn(0);
844   }
845 #endif
846   if ((!s->solution_only && stepnum == 0) || stepnum == s->total_steps) {
847     ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr);
848     ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr);
849     PetscFunctionReturn(0);
850   }
851 
852   if (s->solution_only) {
853 #ifdef PETSC_HAVE_REVOLVE
854     s->rctx->capo = stepnum;
855     shift = 0;
856     whattodo = revolve_action(&s->rctx->check,&s->rctx->capo,&s->rctx->fine,s->rctx->snaps_in,&s->rctx->info,&s->rctx->where);
857     printwhattodo(whattodo,s->rctx,shift);
858 #endif
859   }
860   /* restore a checkpoint */
861   ierr = StackTop(s,&e);CHKERRQ(ierr);
862   if (e && e->stepnum >= stepnum) {
863     SETERRQ2(s->comm,PETSC_ERR_ARG_OUTOFRANGE,"The current step no. is %D, but the step number at top of the stack is %D",stepnum,e->stepnum);
864   }
865   ierr = UpdateTS(ts,s,e);CHKERRQ(ierr);
866   if (s->solution_only || (!s->solution_only && e->stepnum<stepnum)) { /* must recompute */
867     steps = ts->steps;
868     s->recompute = PETSC_TRUE;
869     ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr);
870     ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr);
871     ierr = ReCompute(ts,s,e,stepnum);CHKERRQ(ierr);
872     ts->steps = steps;
873     //ts->total_steps = stepnum;
874     ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr);
875     ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr);
876     #ifdef PETSC_HAVE_REVOLVE
877       s->rctx->reverseonestep = PETSC_FALSE;
878     #endif
879   }
880   PetscFunctionReturn(0);
881 }
882 
883 #undef __FUNCT__
884 #define __FUNCT__ "SetTrajTLR"
885 static PetscErrorCode SetTrajTLR(TS ts,Stack *s,PetscInt stepnum,PetscReal time,Vec X)
886 {
887   PetscInt       whattodo,localstepnum,id;
888   StackElement   e;
889   RevolveCTX     *rctx = s->rctx;
890   PetscErrorCode ierr;
891 
892   PetscFunctionBegin;
893   localstepnum = stepnum%s->stride;
894   if (localstepnum == 0 && stepnum != s->total_steps && !s->recompute) { /* save to disk */
895     id = stepnum/s->stride;
896     if (s->save_stack) {
897       if (stepnum) { /* skip step 0 */
898 #ifdef PETSC_HAVE_REVOLVE
899         PetscPrintf(PETSC_COMM_WORLD,"\x1B[33mDump stack to file\033[0m\n");
900 #endif
901         ierr = StackDumpAll(ts,s,id);CHKERRQ(ierr);
902         s->top = -1; /* reset top */
903 #ifdef PETSC_HAVE_REVOLVE
904         revolve_reset();
905         revolve_create_offline(s->stride,s->max_cps_ram);
906         rctx = s->rctx;
907         rctx->check = 0;
908         rctx->capo  = 0;
909         rctx->fine  = s->stride;
910       }
911 #endif
912     } else {
913       ierr = DumpSingle(ts,s,id);CHKERRQ(ierr);
914       PetscPrintf(PETSC_COMM_WORLD,"\x1B[33mDump a single point to file\033[0m\n");
915     }
916   }
917   /* first forward sweep only checkpoints once in each stride */
918   if (!s->recompute && !s->save_stack) PetscFunctionReturn(0);
919 
920   ierr = ApplyRevolve(ts,s,stepnum,localstepnum,&whattodo);CHKERRQ(ierr);
921   if (whattodo == 2){
922     if (localstepnum < s->top) SETERRQ(s->comm,PETSC_ERR_MEMC,"Illegal modification of a non-top stack element");
923     ierr = ElementCreate(ts,s,&e,localstepnum,time,X);CHKERRQ(ierr);
924     ierr = StackPush(s,e);CHKERRQ(ierr);
925   }
926   PetscFunctionReturn(0);
927 }
928 
929 #undef __FUNCT__
930 #define __FUNCT__ "GetTrajTLR"
931 static PetscErrorCode GetTrajTLR(TS ts,Stack *s,PetscInt stepnum)
932 {
933 #ifdef PETSC_HAVE_REVOLVE
934   PetscInt       whattodo,shift;
935 #endif
936   PetscInt       i,steps,localstepnum,id;
937   PetscReal      stepsize;
938   StackElement   e;
939   PetscErrorCode ierr;
940 
941   PetscFunctionBegin;
942   localstepnum = stepnum%s->stride;
943   if (localstepnum == 0 && stepnum != s->total_steps) { /* load from disk */
944 #ifdef PETSC_HAVE_REVOLVE
945     PetscPrintf(PETSC_COMM_WORLD,"\x1B[33mLoad stack from file\033[0m\n");
946 #endif
947     if (s->save_stack) {
948       id = stepnum/s->stride;
949       ierr = StackLoadAll(ts,s,id);CHKERRQ(ierr);
950       s->top = s->stacksize-1;
951     } else {
952       id = stepnum/s->stride-1;
953       ierr = LoadSingle(ts,s,id);CHKERRQ(ierr);
954     }
955 #ifdef PETSC_HAVE_REVOLVE
956     revolve_reset();
957     revolve_create_offline(s->stride,s->max_cps_ram);
958     s->rctx->check = 0;
959     s->rctx->capo  = 0;
960     s->rctx->fine  = s->stride;
961 #endif
962     if (s->save_stack) {
963 #ifdef PETSC_HAVE_REVOLVE
964       whattodo = 0;
965       while(whattodo!=3) { /* stupid revolve */
966         whattodo = revolve_action(&s->rctx->check,&s->rctx->capo,&s->rctx->fine,s->rctx->snaps_in,&s->rctx->info,&s->rctx->where);
967       }
968 #endif
969     } else {
970       /* save ts context */
971       steps = ts->steps;
972       ts->steps = ts->total_steps;
973       s->recompute = PETSC_TRUE;
974       for (i=ts->steps;i<stepnum;i++) { /* assume fixed step size */
975         ierr = TSTrajectorySet(ts->trajectory,ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
976         ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
977         ierr = TSStep(ts);CHKERRQ(ierr);
978         if (ts->event) {
979           ierr = TSEventMonitor(ts);CHKERRQ(ierr);
980         }
981         if (!ts->steprollback) {
982           ierr = TSPostStep(ts);CHKERRQ(ierr);
983         }
984       }
985       ts->steps = steps;
986       ts->total_steps = stepnum;
987     }
988   }
989 #ifdef PETSC_HAVE_REVOLVE
990   if ( s->rctx->reverseonestep) { /* ready for the reverse step */
991     ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr);
992     ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr);
993     s->rctx->reverseonestep = PETSC_FALSE;
994     PetscFunctionReturn(0);
995   }
996 #endif
997   if ((!s->solution_only && stepnum == 0) || stepnum == s->total_steps) {
998     ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr);
999     ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr);
1000     PetscFunctionReturn(0);
1001   }
1002 
1003   if (s->solution_only) {
1004 #ifdef PETSC_HAVE_REVOLVE
1005     s->rctx->capo = stepnum;
1006     shift = 0;
1007     whattodo = revolve_action(&s->rctx->check,&s->rctx->capo,&s->rctx->fine,s->rctx->snaps_in,&s->rctx->info,&s->rctx->where);
1008     printwhattodo(whattodo,s->rctx,shift);
1009 #endif
1010   }
1011   /* restore a checkpoint */
1012   ierr = StackTop(s,&e);CHKERRQ(ierr);
1013   if (e && e->stepnum >= stepnum) {
1014     SETERRQ2(s->comm,PETSC_ERR_ARG_OUTOFRANGE,"The current step no. is %D, but the step number at top of the stack is %D",stepnum,e->stepnum);
1015   }
1016   ierr = UpdateTS(ts,s,e);CHKERRQ(ierr);
1017   if (s->solution_only || (!s->solution_only && e->stepnum<stepnum)) { /* must recompute */
1018     steps = ts->steps;
1019     s->recompute = PETSC_TRUE;
1020     ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr);
1021     ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr);
1022     ierr = ReCompute(ts,s,e,stepnum);CHKERRQ(ierr);
1023     ts->steps = steps;
1024     //ts->total_steps = stepnum;
1025     ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr);
1026     ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr);
1027     #ifdef PETSC_HAVE_REVOLVE
1028       s->rctx->reverseonestep = PETSC_FALSE;
1029     #endif
1030   } else {
1031     ierr = ElementDestroy(s,e);CHKERRQ(ierr);
1032   }
1033   PetscFunctionReturn(0);
1034 }
1035 
1036 #undef __FUNCT__
1037 #define __FUNCT__ "SetTrajTLNR"
1038 static PetscErrorCode SetTrajTLNR(TS ts,Stack *s,PetscInt stepnum,PetscReal time,Vec X)
1039 {
1040   PetscInt       whattodo,localstepnum,id;
1041   StackElement   e;
1042   PetscErrorCode ierr;
1043 
1044   PetscFunctionBegin;
1045   localstepnum = stepnum%s->stride;
1046   if (localstepnum == 0 && stepnum != s->total_steps && stepnum != 0  && !s->recompute) {
1047     id = stepnum/s->stride;
1048     if (s->save_stack) {
1049       ierr = StackDumpAll(ts,s,id);CHKERRQ(ierr);
1050       s->top = -1; /* reset top */
1051     } else {
1052       ierr = DumpSingle(ts,s,id);CHKERRQ(ierr);
1053     }
1054   }
1055 
1056   ierr = ApplyRevolve(ts,s,stepnum,localstepnum,&whattodo);CHKERRQ(ierr);
1057   if (whattodo == 2){
1058     if (localstepnum < s->top) SETERRQ(s->comm,PETSC_ERR_MEMC,"Illegal modification of a non-top stack element");
1059     ierr = ElementCreate(ts,s,&e,localstepnum,time,X);CHKERRQ(ierr);
1060     ierr = StackPush(s,e);CHKERRQ(ierr);
1061   }
1062   PetscFunctionReturn(0);
1063 }
1064 
1065 #undef __FUNCT__
1066 #define __FUNCT__ "GetTrajTLNR"
1067 static PetscErrorCode GetTrajTLNR(TS ts,Stack *s,PetscInt stepnum)
1068 {
1069   PetscInt       steps,id,localstepnum;
1070   PetscReal      stepsize;
1071   StackElement   e;
1072   PetscErrorCode ierr;
1073 
1074   PetscFunctionBegin;
1075   localstepnum = stepnum%s->stride;
1076   if (stepnum != s->total_steps && localstepnum==0) {
1077     id = stepnum/s->stride;
1078     if (s->save_stack) {
1079       ierr = StackLoadAll(ts,s,id);CHKERRQ(ierr);
1080     } else {
1081       ierr = LoadSingle(ts,s,id);CHKERRQ(ierr);
1082     }
1083   }
1084   if ((!s->solution_only && stepnum == 0) || stepnum == s->total_steps) {
1085     ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr);
1086     ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr);
1087     PetscFunctionReturn(0);
1088   }
1089 
1090   /* restore a checkpoint */
1091   ierr = StackTop(s,&e);CHKERRQ(ierr);
1092   if (e && e->stepnum >= stepnum) {
1093     SETERRQ2(s->comm,PETSC_ERR_ARG_OUTOFRANGE,"The current step no. is %D, but the step number at top of the stack is %D",stepnum,e->stepnum);
1094   }
1095   ierr = UpdateTS(ts,s,e);CHKERRQ(ierr);
1096   if (s->solution_only || (!s->solution_only && e->stepnum<stepnum)) { /* must recompute */
1097     steps = ts->steps;
1098     s->recompute = PETSC_TRUE;
1099     ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr);
1100     ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr);
1101     ierr = ReCompute(ts,s,e,stepnum);CHKERRQ(ierr);
1102     ts->steps = steps;
1103     //ts->total_steps = stepnum;
1104     ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr);
1105     ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr);
1106   } else {
1107     ierr = ElementDestroy(s,e);CHKERRQ(ierr);
1108   }
1109   PetscFunctionReturn(0);
1110 }
1111 
1112 
1113 #undef __FUNCT__
1114 #define __FUNCT__ "SetTrajRMS"
1115 static PetscErrorCode SetTrajRMS(TS ts,Stack *s,PetscInt stepnum,PetscReal time,Vec X)
1116 {
1117   PetscInt       whattodo;
1118   StackElement   e;
1119   PetscErrorCode ierr;
1120 
1121   PetscFunctionBegin;
1122   ierr = ApplyRevolve(ts,s,stepnum,stepnum,&whattodo);CHKERRQ(ierr);
1123   if (whattodo == 2){
1124     if (stepnum < s->top) SETERRQ(s->comm,PETSC_ERR_MEMC,"Illegal modification of a non-top stack element");
1125     ierr = ElementCreate(ts,s,&e,stepnum,time,X);CHKERRQ(ierr);
1126     ierr = StackPush(s,e);CHKERRQ(ierr);
1127   }
1128   PetscFunctionReturn(0);
1129 }
1130 
1131 #undef __FUNCT__
1132 #define __FUNCT__ "GetTrajRMS"
1133 static PetscErrorCode GetTrajRMS(TS ts,Stack *s,PetscInt stepnum)
1134 {
1135 #ifdef PETSC_HAVE_REVOLVE
1136   PetscInt       whattodo,shift;
1137 #endif
1138   PetscInt       steps;
1139   PetscReal      stepsize;
1140   StackElement   e;
1141   PetscErrorCode ierr;
1142 
1143   PetscFunctionBegin;
1144 #ifdef PETSC_HAVE_REVOLVE
1145   if ( s->rctx->reverseonestep) { /* ready for the reverse step */
1146     ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr);
1147     ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr);
1148     s->rctx->reverseonestep = PETSC_FALSE;
1149     PetscFunctionReturn(0);
1150   }
1151 #endif
1152   if ((!s->solution_only && stepnum == 0) || stepnum == s->total_steps) {
1153     ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr);
1154     ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr);
1155     PetscFunctionReturn(0);
1156   }
1157 
1158   if (s->solution_only) {
1159 #ifdef PETSC_HAVE_REVOLVE
1160     s->rctx->capo = stepnum;
1161     shift = 0;
1162     whattodo = revolve_action(&s->rctx->check,&s->rctx->capo,&s->rctx->fine,s->rctx->snaps_in,&s->rctx->info,&s->rctx->where);
1163     printwhattodo(whattodo,s->rctx,shift);
1164 #endif
1165   }
1166   /* restore a checkpoint */
1167   if (!s->rctx->where) {
1168     ierr = LoadSingle(ts,s,stepnum);CHKERRQ(ierr);
1169   } else {
1170     ierr = StackTop(s,&e);CHKERRQ(ierr);
1171     if (e && e->stepnum >= stepnum) {
1172       SETERRQ2(s->comm,PETSC_ERR_ARG_OUTOFRANGE,"The current step no. is %D, but the step number at top of the stack is %D",stepnum,e->stepnum);
1173     }
1174     ierr = UpdateTS(ts,s,e);CHKERRQ(ierr);
1175   }
1176 
1177   if (s->solution_only || (!s->solution_only && e->stepnum<stepnum)) { /* must recompute */
1178     steps = ts->steps;
1179     s->recompute = PETSC_TRUE;
1180     ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr);
1181     ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr);
1182     ierr = ReCompute(ts,s,e,stepnum);CHKERRQ(ierr);
1183     ts->steps = steps;
1184     //ts->total_steps = stepnum;
1185     ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr);
1186     ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr);
1187     #ifdef PETSC_HAVE_REVOLVE
1188       s->rctx->reverseonestep = PETSC_FALSE;
1189     #endif
1190   } else {
1191     ierr = ElementDestroy(s,e);CHKERRQ(ierr);
1192   }
1193   PetscFunctionReturn(0);
1194 }
1195 
1196 #undef __FUNCT__
1197 #define __FUNCT__ "TSTrajectorySet_Memory"
1198 PetscErrorCode TSTrajectorySet_Memory(TSTrajectory tj,TS ts,PetscInt stepnum,PetscReal time,Vec X)
1199 {
1200   Stack          *s = (Stack*)tj->data;
1201   PetscErrorCode ierr;
1202 
1203   PetscFunctionBegin;
1204   if (!s->recompute) { /* use global stepnum in the forward sweep */
1205     ierr = TSGetTotalSteps(ts,&stepnum);CHKERRQ(ierr);
1206   }
1207   /* for consistency */
1208   if (!s->recompute && stepnum == 0) ts->ptime_prev = ts->ptime-ts->time_step;
1209   switch (s->stype) {
1210     case NONE:
1211       ierr = SetTrajN(ts,s,stepnum,time,X);CHKERRQ(ierr);
1212       break;
1213     case TWO_LEVEL_NOREVOLVE:
1214       ierr = SetTrajTLNR(ts,s,stepnum,time,X);CHKERRQ(ierr);
1215       break;
1216     case TWO_LEVEL_REVOLVE:
1217       ierr = SetTrajTLR(ts,s,stepnum,time,X);CHKERRQ(ierr);
1218       break;
1219     case REVOLVE_OFFLINE:
1220       ierr = SetTrajROF(ts,s,stepnum,time,X);CHKERRQ(ierr);
1221       break;
1222     case REVOLVE_ONLINE:
1223       ierr = SetTrajRON(ts,s,stepnum,time,X);CHKERRQ(ierr);
1224       break;
1225     case REVOLVE_MULTISTAGE:
1226       ierr = SetTrajRMS(ts,s,stepnum,time,X);CHKERRQ(ierr);
1227       break;
1228     default:
1229       break;
1230   }
1231   PetscFunctionReturn(0);
1232 }
1233 
1234 #undef __FUNCT__
1235 #define __FUNCT__ "TSTrajectoryGet_Memory"
1236 PetscErrorCode TSTrajectoryGet_Memory(TSTrajectory tj,TS ts,PetscInt stepnum,PetscReal *t)
1237 {
1238   Stack          *s = (Stack*)tj->data;
1239   PetscErrorCode ierr;
1240 
1241   PetscFunctionBegin;
1242   ierr = TSGetTotalSteps(ts,&stepnum);CHKERRQ(ierr);
1243   if (stepnum == 0) PetscFunctionReturn(0);
1244   switch (s->stype) {
1245     case NONE:
1246       ierr = GetTrajN(ts,s,stepnum);CHKERRQ(ierr);
1247       break;
1248     case TWO_LEVEL_NOREVOLVE:
1249       ierr = GetTrajTLNR(ts,s,stepnum);CHKERRQ(ierr);
1250       break;
1251     case TWO_LEVEL_REVOLVE:
1252       ierr = GetTrajTLR(ts,s,stepnum);CHKERRQ(ierr);
1253       break;
1254     case REVOLVE_OFFLINE:
1255       ierr = GetTrajROF(ts,s,stepnum);CHKERRQ(ierr);
1256       break;
1257     case REVOLVE_ONLINE:
1258       ierr = GetTrajRON(ts,s,stepnum);CHKERRQ(ierr);
1259       break;
1260     case REVOLVE_MULTISTAGE:
1261       ierr = GetTrajRMS(ts,s,stepnum);CHKERRQ(ierr);
1262       break;
1263     default:
1264       break;
1265   }
1266   PetscFunctionReturn(0);
1267 }
1268 
1269 #undef __FUNCT__
1270 #define __FUNCT__ "TSTrajectoryDestroy_Memory"
1271 PETSC_EXTERN PetscErrorCode TSTrajectoryDestroy_Memory(TSTrajectory tj)
1272 {
1273   Stack          *s = (Stack*)tj->data;
1274   PetscErrorCode ierr;
1275 
1276   PetscFunctionBegin;
1277   if (s->stype > TWO_LEVEL_NOREVOLVE) {
1278 #ifdef PETSC_HAVE_REVOLVE
1279     revolve_reset();
1280 #endif
1281   }
1282   ierr = StackDestroy(s);CHKERRQ(ierr);
1283   PetscFunctionReturn(0);
1284 }
1285 
1286 /*MC
1287       TSTRAJECTORYMEMORY - Stores each solution of the ODE/ADE in memory
1288 
1289   Level: intermediate
1290 
1291 .seealso:  TSTrajectoryCreate(), TS, TSTrajectorySetType()
1292 
1293 M*/
1294 #undef __FUNCT__
1295 #define __FUNCT__ "TSTrajectoryCreate_Memory"
1296 PETSC_EXTERN PetscErrorCode TSTrajectoryCreate_Memory(TSTrajectory tj,TS ts)
1297 {
1298   Stack          *s;
1299   PetscErrorCode ierr;
1300 
1301   PetscFunctionBegin;
1302   tj->ops->set            = TSTrajectorySet_Memory;
1303   tj->ops->get            = TSTrajectoryGet_Memory;
1304   tj->ops->setup          = TSTrajectorySetUp_Memory;
1305   tj->ops->destroy        = TSTrajectoryDestroy_Memory;
1306   tj->ops->setfromoptions = TSTrajectorySetFromOptions_Memory;
1307 
1308   ierr = PetscCalloc1(1,&s);CHKERRQ(ierr);
1309   s->stype        = NONE;
1310   s->max_cps_ram  = -1; /* -1 indicates that it is not set */
1311   s->max_cps_disk = -1; /* -1 indicates that it is not set */
1312   s->stride       = 0; /* if not zero, two-level checkpointing will be used */
1313   s->use_online   = PETSC_FALSE;
1314   s->save_stack   = PETSC_TRUE;
1315   s->solution_only= PETSC_TRUE;
1316   tj->data        = s;
1317   PetscFunctionReturn(0);
1318 }
1319