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