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