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