xref: /petsc/src/ts/trajectory/impls/memory/trajmemory.c (revision 164d268a552a6708f7fd68c243172cce5abc720b)
1 #include <petsc/private/tsimpl.h>        /*I "petscts.h"  I*/
2 #include <petscsys.h>
3 
4 #ifdef PETSC_HAVE_REVOLVE
5 #include <revolve_c.h>
6 #endif
7 
8 PetscLogEvent Disk_Write, Disk_Read;
9 
10 typedef enum {NONE,TWO_LEVEL_NOREVOLVE,TWO_LEVEL_REVOLVE,REVOLVE_OFFLINE,REVOLVE_ONLINE,REVOLVE_MULTISTAGE} SchedulerType;
11 
12 typedef struct _StackElement {
13   PetscInt  stepnum;
14   Vec       X;
15   Vec       *Y;
16   PetscReal time;
17   PetscReal timeprev;
18 } *StackElement;
19 
20 typedef struct _RevolveCTX {
21   PetscBool reverseonestep;
22   PetscInt  where;
23   PetscInt  snaps_in;
24   PetscInt  stepsleft;
25   PetscInt  check;
26   PetscInt  oldcapo;
27   PetscInt  capo;
28   PetscInt  fine;
29   PetscInt  info;
30 } RevolveCTX;
31 
32 typedef struct _Stack {
33   SchedulerType stype;
34   PetscBool     use_online;
35   PetscBool     recompute;
36   PetscBool     solution_only;
37   PetscBool     save_stack;
38   MPI_Comm      comm;
39   RevolveCTX    *rctx;
40   PetscInt      max_cps_ram;  /* maximum checkpoints in RAM */
41   PetscInt      max_cps_disk; /* maximum checkpoints on disk */
42   PetscInt      stride;
43   PetscInt      total_steps;  /* total number of steps */
44   PetscInt      numY;
45   PetscInt      stacksize;
46   PetscInt      top;          /* top of the stack */
47   StackElement  *stack;       /* container */
48 } Stack;
49 
50 static PetscErrorCode StackCreate(MPI_Comm,Stack *,PetscInt,PetscInt);
51 static PetscErrorCode StackDestroy(Stack*);
52 static PetscErrorCode StackPush(Stack*,StackElement);
53 static PetscErrorCode StackPop(Stack*,StackElement*);
54 static PetscErrorCode StackTop(Stack*,StackElement*);
55 static PetscErrorCode StackDumpAll(TS,Stack*,PetscInt);
56 static PetscErrorCode StackLoadAll(TS,Stack*,PetscInt);
57 
58 #ifdef PETSC_HAVE_REVOLVE
59 static void printwhattodo(PetscInt whattodo,RevolveCTX *rctx,PetscInt shift)
60 {
61   switch(whattodo) {
62     case 1:
63       PetscPrintf(PETSC_COMM_WORLD,"\x1B[35mAdvance from %D to %D\033[0m\n",rctx->oldcapo+shift,rctx->capo+shift);
64       break;
65     case 2:
66       PetscPrintf(PETSC_COMM_WORLD,"\x1B[35mStore in checkpoint number %D (located in RAM)\033[0m\n",rctx->check);
67       break;
68     case 3:
69       PetscPrintf(PETSC_COMM_WORLD,"\x1B[35mFirst turn: Initialize adjoints and reverse first step\033[0m\n");
70       break;
71     case 4:
72       PetscPrintf(PETSC_COMM_WORLD,"\x1B[35mForward and reverse one step\033[0m\n");
73       break;
74     case 5:
75       PetscPrintf(PETSC_COMM_WORLD,"\x1B[35mRestore in checkpoint number %D (located in RAM)\033[0m\n",rctx->check);
76       break;
77     case 7:
78       PetscPrintf(PETSC_COMM_WORLD,"\x1B[35mStore in checkpoint number %D (located on disk)\033[0m\n",rctx->check);
79       break;
80     case 8:
81       PetscPrintf(PETSC_COMM_WORLD,"\x1B[35mRestore in checkpoint number %D (located on disk)\033[0m\n",rctx->check);
82       break;
83     case -1:
84       PetscPrintf(PETSC_COMM_WORLD,"\x1B[35mError!");
85       break;
86   }
87 }
88 #endif
89 
90 #undef __FUNCT__
91 #define __FUNCT__ "StackCreate"
92 static PetscErrorCode StackCreate(MPI_Comm comm,Stack *s,PetscInt size,PetscInt ny)
93 {
94   PetscErrorCode ierr;
95 
96   PetscFunctionBegin;
97   s->top  = -1;
98   s->comm = comm;
99   s->numY = ny;
100 
101   ierr = PetscMalloc1(size*sizeof(StackElement),&s->stack);CHKERRQ(ierr);
102   PetscFunctionReturn(0);
103 }
104 
105 #undef __FUNCT__
106 #define __FUNCT__ "StackDestroy"
107 static PetscErrorCode StackDestroy(Stack *s)
108 {
109   PetscInt       i;
110   PetscErrorCode ierr;
111 
112   PetscFunctionBegin;
113   if (s->top>-1) {
114     for (i=0;i<=s->top;i++) {
115       ierr = VecDestroy(&s->stack[i]->X);CHKERRQ(ierr);
116       if (!s->solution_only) {
117         ierr = VecDestroyVecs(s->numY,&s->stack[i]->Y);CHKERRQ(ierr);
118       }
119       ierr = PetscFree(s->stack[i]);CHKERRQ(ierr);
120     }
121   }
122   ierr = PetscFree(s->stack);CHKERRQ(ierr);
123   if (s->stype) {
124     ierr = PetscFree(s->rctx);CHKERRQ(ierr);
125   }
126   ierr = PetscFree(s);CHKERRQ(ierr);
127   PetscFunctionReturn(0);
128 }
129 
130 #undef __FUNCT__
131 #define __FUNCT__ "StackPush"
132 static PetscErrorCode StackPush(Stack *s,StackElement e)
133 {
134   PetscFunctionBegin;
135   if (s->top+1 >= s->stacksize) SETERRQ1(s->comm,PETSC_ERR_MEMC,"Maximum stack size (%D) exceeded",s->stacksize);
136   s->stack[++s->top] = e;
137   PetscFunctionReturn(0);
138 }
139 
140 #undef __FUNCT__
141 #define __FUNCT__ "StackPop"
142 static PetscErrorCode StackPop(Stack *s,StackElement *e)
143 {
144   PetscFunctionBegin;
145   if (s->top == -1) SETERRQ(s->comm,PETSC_ERR_MEMC,"Empty stack");
146   *e = s->stack[s->top--];
147   PetscFunctionReturn(0);
148 }
149 
150 #undef __FUNCT__
151 #define __FUNCT__ "StackTop"
152 static PetscErrorCode StackTop(Stack *s,StackElement *e)
153 {
154   PetscFunctionBegin;
155   *e = s->stack[s->top];
156   PetscFunctionReturn(0);
157 }
158 
159 #undef __FUNCT__
160 #define __FUNCT__ "StackFind"
161 static PetscErrorCode StackFind(Stack *s,StackElement *e,PetscInt index)
162 {
163   PetscFunctionBegin;
164   *e = s->stack[index];
165   PetscFunctionReturn(0);
166 }
167 
168 #undef __FUNCT__
169 #define __FUNCT__ "OutputBIN"
170 static PetscErrorCode OutputBIN(const char *filename,PetscViewer *viewer)
171 {
172   PetscErrorCode ierr;
173 
174   PetscFunctionBegin;
175   ierr = PetscViewerCreate(PETSC_COMM_WORLD,viewer);CHKERRQ(ierr);
176   ierr = PetscViewerSetType(*viewer,PETSCVIEWERBINARY);CHKERRQ(ierr);
177   ierr = PetscViewerFileSetMode(*viewer,FILE_MODE_WRITE);CHKERRQ(ierr);
178   ierr = PetscViewerFileSetName(*viewer,filename);CHKERRQ(ierr);
179   PetscFunctionReturn(0);
180 }
181 
182 #undef __FUNCT__
183 #define __FUNCT__ "WriteToDisk"
184 static PetscErrorCode WriteToDisk(PetscInt stepnum,PetscReal time,PetscReal timeprev,Vec X,Vec *Y,PetscInt numY,PetscBool solution_only,PetscViewer viewer)
185 {
186   PetscInt       i;
187   PetscErrorCode ierr;
188 
189   PetscFunctionBegin;
190   ierr = PetscViewerBinaryWrite(viewer,&stepnum,1,PETSC_INT,PETSC_FALSE);CHKERRQ(ierr);
191   ierr = VecView(X,viewer);CHKERRQ(ierr);
192   for (i=0;!solution_only && i<numY;i++) {
193     ierr = VecView(Y[i],viewer);CHKERRQ(ierr);
194   }
195   ierr = PetscViewerBinaryWrite(viewer,&time,1,PETSC_REAL,PETSC_FALSE);CHKERRQ(ierr);
196   ierr = PetscViewerBinaryWrite(viewer,&timeprev,1,PETSC_REAL,PETSC_FALSE);CHKERRQ(ierr);
197   PetscFunctionReturn(0);
198 }
199 
200 #undef __FUNCT__
201 #define __FUNCT__ "ReadFromDisk"
202 static PetscErrorCode ReadFromDisk(PetscInt *stepnum,PetscReal *time,PetscReal *timeprev,Vec X,Vec *Y,PetscInt numY,PetscBool solution_only,PetscViewer viewer)
203 {
204   PetscInt       i;
205   PetscErrorCode ierr;
206 
207   PetscFunctionBegin;
208   ierr = PetscViewerBinaryRead(viewer,stepnum,1,NULL,PETSC_INT);CHKERRQ(ierr);
209   ierr = VecLoad(X,viewer);CHKERRQ(ierr);
210   for (i=0;!solution_only && i<numY;i++) {
211     ierr = VecLoad(Y[i],viewer);CHKERRQ(ierr);
212   }
213   ierr = PetscViewerBinaryRead(viewer,time,1,NULL,PETSC_REAL);CHKERRQ(ierr);
214   ierr = PetscViewerBinaryRead(viewer,timeprev,1,NULL,PETSC_REAL);CHKERRQ(ierr);
215   PetscFunctionReturn(0);
216 }
217 
218 #undef __FUNCT__
219 #define __FUNCT__ "StackDumpAll"
220 static PetscErrorCode StackDumpAll(TS ts,Stack *s,PetscInt id)
221 {
222   Vec            *Y;
223   PetscInt       i;
224   StackElement   e;
225   PetscViewer    viewer;
226   char           filename[PETSC_MAX_PATH_LEN];
227   PetscErrorCode ierr;
228 
229   PetscFunctionBegin;
230   if (id == 1) {
231     PetscMPIInt rank;
232     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)ts),&rank);CHKERRQ(ierr);
233     if (!rank) {
234       ierr = PetscRMTree("SA-data");CHKERRQ(ierr);
235       ierr = PetscMkdir("SA-data");CHKERRQ(ierr);
236     }
237   }
238   ierr = PetscSNPrintf(filename,sizeof(filename),"SA-data/SA-STACK%06d.bin",id);CHKERRQ(ierr);
239   ierr = OutputBIN(filename,&viewer);CHKERRQ(ierr);
240   for (i=0;i<s->stacksize;i++) {
241     e = s->stack[i];
242     ierr = WriteToDisk(e->stepnum,e->time,e->timeprev,e->X,e->Y,s->numY,s->solution_only,viewer);CHKERRQ(ierr);
243   }
244   /* save the last step for restart, the last step is in memory when using single level schemes, but not necessarily the case for multi level schemes */
245   ierr = TSGetStages(ts,&s->numY,&Y);CHKERRQ(ierr);
246   ierr = WriteToDisk(ts->total_steps,ts->ptime,ts->ptime_prev,ts->vec_sol,Y,s->numY,s->solution_only,viewer);CHKERRQ(ierr);
247   for (i=0;i<s->stacksize;i++) {
248     ierr = StackPop(s,&e);CHKERRQ(ierr);
249     ierr = VecDestroy(&e->X);CHKERRQ(ierr);
250     if (!s->solution_only) {
251       ierr = VecDestroyVecs(s->numY,&e->Y);CHKERRQ(ierr);
252     }
253     ierr = PetscFree(e);CHKERRQ(ierr);
254   }
255   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
256   PetscFunctionReturn(0);
257 }
258 
259 #undef __FUNCT__
260 #define __FUNCT__ "StackLoadAll"
261 static PetscErrorCode StackLoadAll(TS ts,Stack *s,PetscInt id)
262 {
263   Vec            *Y;
264   PetscInt       i;
265   StackElement   e;
266   PetscViewer    viewer;
267   char           filename[PETSC_MAX_PATH_LEN];
268   PetscErrorCode ierr;
269 
270   PetscFunctionBegin;
271   ierr = PetscSNPrintf(filename,sizeof filename,"SA-data/SA-STACK%06d.bin",id);CHKERRQ(ierr);
272   ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);CHKERRQ(ierr);
273   for (i=0;i<s->stacksize;i++) {
274     ierr = PetscCalloc1(1,&e);
275     ierr = TSGetStages(ts,&s->numY,&Y);CHKERRQ(ierr);
276     ierr = VecDuplicate(Y[0],&e->X);CHKERRQ(ierr);
277     if (!s->solution_only && s->numY>0) {
278       ierr = VecDuplicateVecs(Y[0],s->numY,&e->Y);CHKERRQ(ierr);
279     }
280     ierr = StackPush(s,e);CHKERRQ(ierr);
281   }
282   for (i=0;i<s->stacksize;i++) {
283     e = s->stack[i];
284     ierr = ReadFromDisk(&e->stepnum,&e->time,&e->timeprev,e->X,e->Y,s->numY,s->solution_only,viewer);CHKERRQ(ierr);
285   }
286   /* load the last step into TS */
287   ierr = TSGetStages(ts,&s->numY,&Y);CHKERRQ(ierr);
288   ierr = ReadFromDisk(&ts->total_steps,&ts->ptime,&ts->ptime_prev,ts->vec_sol,Y,s->numY,s->solution_only,viewer);CHKERRQ(ierr);
289   ierr = TSSetTimeStep(ts,ts->ptime-ts->ptime_prev);CHKERRQ(ierr);
290   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
291   PetscFunctionReturn(0);
292 }
293 
294 #undef __FUNCT__
295 #define __FUNCT__ "DumpSingle"
296 static PetscErrorCode DumpSingle(TS ts,Stack *s,PetscInt id)
297 {
298   Vec            *Y;
299   PetscInt       stepnum;
300   PetscViewer    viewer;
301   char           filename[PETSC_MAX_PATH_LEN];
302   PetscErrorCode ierr;
303 
304   PetscFunctionBegin;
305   ierr = TSGetTotalSteps(ts,&stepnum);CHKERRQ(ierr);
306   if (id == 0) {
307     PetscMPIInt rank;
308     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)ts),&rank);CHKERRQ(ierr);
309     if (!rank) {
310       ierr = PetscRMTree("SA-data");CHKERRQ(ierr);
311       ierr = PetscMkdir("SA-data");CHKERRQ(ierr);
312     }
313   }
314   ierr = PetscSNPrintf(filename,sizeof(filename),"SA-data/SA-CPS%06d.bin",id);CHKERRQ(ierr);
315   ierr = OutputBIN(filename,&viewer);CHKERRQ(ierr);
316 
317   ierr = PetscLogEventBegin(Disk_Write,ts,0,0,0);CHKERRQ(ierr);
318   ierr = TSGetStages(ts,&s->numY,&Y);CHKERRQ(ierr);
319   ierr = WriteToDisk(stepnum,ts->ptime,ts->ptime_prev,ts->vec_sol,Y,s->numY,s->solution_only,viewer);CHKERRQ(ierr);
320   ierr = PetscLogEventEnd(Disk_Write,ts,0,0,0);CHKERRQ(ierr);
321 
322   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
323   PetscFunctionReturn(0);
324 }
325 
326 #undef __FUNCT__
327 #define __FUNCT__ "LoadSingle"
328 static PetscErrorCode LoadSingle(TS ts,Stack *s,PetscInt id)
329 {
330   Vec            *Y;
331   PetscViewer    viewer;
332   char           filename[PETSC_MAX_PATH_LEN];
333   PetscErrorCode ierr;
334 
335   PetscFunctionBegin;
336   ierr = PetscSNPrintf(filename,sizeof filename,"SA-data/SA-CPS%06d.bin",id);CHKERRQ(ierr);
337   ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);CHKERRQ(ierr);
338 
339   ierr = PetscLogEventBegin(Disk_Read,ts,0,0,0);CHKERRQ(ierr);
340   ierr = TSGetStages(ts,&s->numY,&Y);CHKERRQ(ierr);
341   ierr = ReadFromDisk(&ts->total_steps,&ts->ptime,&ts->ptime_prev,ts->vec_sol,Y,s->numY,s->solution_only,viewer);CHKERRQ(ierr);
342   ierr = PetscLogEventEnd(Disk_Read,ts,0,0,0);CHKERRQ(ierr);
343 
344   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
345   PetscFunctionReturn(0);
346 }
347 
348 #undef __FUNCT__
349 #define __FUNCT__ "TSTrajectorySetStride_Memory"
350 PetscErrorCode TSTrajectorySetStride_Memory(TSTrajectory tj,TS ts,PetscInt stride)
351 {
352   Stack *s = (Stack*)tj->data;
353 
354   PetscFunctionBegin;
355   s->stride = stride;
356   PetscFunctionReturn(0);
357 }
358 
359 #undef __FUNCT__
360 #define __FUNCT__ "TSTrajectorySetMaxCpsRAM_Memory"
361 PetscErrorCode TSTrajectorySetMaxCpsRAM_Memory(TSTrajectory tj,TS ts,PetscInt max_cps_ram)
362 {
363   Stack *s = (Stack*)tj->data;
364 
365   PetscFunctionBegin;
366   s->max_cps_ram = max_cps_ram;
367   PetscFunctionReturn(0);
368 }
369 
370 #undef __FUNCT__
371 #define __FUNCT__ "TSTrajectorySetMaxCpsDisk_Memory"
372 PetscErrorCode TSTrajectorySetMaxCpsDisk_Memory(TSTrajectory tj,TS ts,PetscInt max_cps_disk)
373 {
374   Stack *s = (Stack*)tj->data;
375 
376   PetscFunctionBegin;
377   s->max_cps_disk = max_cps_disk;
378   PetscFunctionReturn(0);
379 }
380 
381 #undef __FUNCT__
382 #define __FUNCT__ "TSTrajectorySetRevolveOnline"
383 PetscErrorCode TSTrajectorySetRevolveOnline(TSTrajectory tj,PetscBool use_online)
384 {
385   Stack *s = (Stack*)tj->data;
386 
387   PetscFunctionBegin;
388   s->use_online = use_online;
389   PetscFunctionReturn(0);
390 }
391 
392 #undef __FUNCT__
393 #define __FUNCT__ "TSTrajectorySetSaveStack"
394 PetscErrorCode TSTrajectorySetSaveStack(TSTrajectory tj,PetscBool save_stack)
395 {
396   Stack *s = (Stack*)tj->data;
397 
398   PetscFunctionBegin;
399   s->save_stack = save_stack;
400   PetscFunctionReturn(0);
401 }
402 
403 #undef __FUNCT__
404 #define __FUNCT__ "TSTrajectorySetSolutionOnly"
405 PetscErrorCode TSTrajectorySetSolutionOnly(TSTrajectory tj,PetscBool solution_only)
406 {
407   Stack *s = (Stack*)tj->data;
408 
409   PetscFunctionBegin;
410   s->solution_only = solution_only;
411   PetscFunctionReturn(0);
412 }
413 
414 #undef __FUNCT__
415 #define __FUNCT__ "TSTrajectorySetFromOptions_Memory"
416 PetscErrorCode TSTrajectorySetFromOptions_Memory(PetscOptions *PetscOptionsObject,TSTrajectory tj)
417 {
418   Stack          *s = (Stack*)tj->data;
419   PetscErrorCode ierr;
420 
421   PetscFunctionBegin;
422   ierr = PetscOptionsHead(PetscOptionsObject,"Memory based TS trajectory options");CHKERRQ(ierr);
423   {
424     ierr = PetscOptionsInt("-tstrajectory_max_cps_ram","Maximum number of checkpoints in RAM","TSTrajectorySetMaxCpsRAM_Memory",s->max_cps_ram,&s->max_cps_ram,NULL);CHKERRQ(ierr);
425     ierr = PetscOptionsInt("-tstrajectory_max_cps_disk","Maximum number of checkpoints on disk","TSTrajectorySetMaxCpsDisk_Memory",s->max_cps_disk,&s->max_cps_disk,NULL);CHKERRQ(ierr);
426     ierr = PetscOptionsInt("-tstrajectory_stride","Stride to save checkpoints to file","TSTrajectorySetStride_Memory",s->stride,&s->stride,NULL);CHKERRQ(ierr);
427     ierr = PetscOptionsBool("-tstrajectory_revolve_online","Trick TS trajectory into using online mode of revolve","TSTrajectorySetRevolveOnline",s->use_online,&s->use_online,NULL);CHKERRQ(ierr);
428     ierr = PetscOptionsBool("-tstrajectory_save_stack","Save all stack to disk","TSTrajectorySetSaveStack",s->save_stack,&s->save_stack,NULL);CHKERRQ(ierr);
429     ierr = PetscOptionsBool("-tstrajectory_solution_only","Checkpoint solution only","TSTrajectorySetSolutionOnly",s->solution_only,&s->solution_only,NULL);CHKERRQ(ierr);
430   }
431   ierr = PetscOptionsTail();CHKERRQ(ierr);
432   PetscFunctionReturn(0);
433 }
434 
435 #undef __FUNCT__
436 #define __FUNCT__ "TSTrajectorySetUp_Memory"
437 PetscErrorCode TSTrajectorySetUp_Memory(TSTrajectory tj,TS ts)
438 {
439   Stack          *s = (Stack*)tj->data;
440   RevolveCTX     *rctx;
441   PetscInt       numY;
442   PetscBool      flg;
443   PetscErrorCode ierr;
444 
445   PetscFunctionBegin;
446   PetscStrcmp(((PetscObject)ts->adapt)->type_name,TSADAPTNONE,&flg);
447   if (flg) { /* fixed time step */
448     s->total_steps = PetscMin(ts->max_steps,(PetscInt)(ceil(ts->max_time/ts->time_step)));
449   }
450   if (s->max_cps_ram > 1) s->stacksize = s->max_cps_ram;
451   if (s->stride > 1) { /* two level mode works for both fixed time step and adaptive time step */
452     if (s->max_cps_ram > 1 && s->max_cps_ram < s->stride-1) { /* use revolve_offline for each stride */
453       s->stype = TWO_LEVEL_REVOLVE;
454     }else { /* checkpoint all for each stride */
455       s->stype     = TWO_LEVEL_NOREVOLVE;
456       s->stacksize = s->solution_only ? s->stride : s->stride-1;
457     }
458   } else {
459     if (flg) { /* fixed time step */
460       if (s->max_cps_ram >= s->total_steps-1 || s->max_cps_ram < 1) { /* checkpoint all */
461         s->stype     = NONE;
462         s->stacksize = s->solution_only ? s->total_steps : s->total_steps-1;
463       } else {
464         if (s->max_cps_disk > 1) { /* disk can be used */
465           s->stype = REVOLVE_MULTISTAGE;
466         } else { /* memory only */
467           s->stype = REVOLVE_OFFLINE;
468         }
469       }
470     } else { /* adaptive time step */
471       s->stype = REVOLVE_ONLINE;
472     }
473     if (s->use_online) { /* trick into online */
474       s->stype     = REVOLVE_ONLINE;
475       s->stacksize = s->max_cps_ram;
476     }
477   }
478 
479   if (s->stype > TWO_LEVEL_NOREVOLVE) {
480 #ifndef PETSC_HAVE_REVOLVE
481     SETERRQ(s->comm,PETSC_ERR_SUP,"revolve is needed when there is not enough memory to checkpoint all time steps according to the user's settings, please reconfigure with the additional option --download-revolve.");
482 #else
483     if (s->stype == TWO_LEVEL_REVOLVE) revolve_create_offline(s->stride,s->max_cps_ram);
484     else if (s->stype == REVOLVE_OFFLINE) revolve_create_offline(s->total_steps,s->max_cps_ram);
485     else if (s->stype == REVOLVE_ONLINE) revolve_create_online(s->max_cps_ram);
486     else if (s->stype ==REVOLVE_MULTISTAGE) revolve_create_multistage(s->total_steps,s->max_cps_ram+s->max_cps_disk,s->max_cps_ram);
487 
488     ierr = PetscCalloc1(1,&rctx);CHKERRQ(ierr);
489     rctx->snaps_in       = s->max_cps_ram; /* for theta methods snaps_in=2*max_cps_ram */
490     rctx->reverseonestep = PETSC_FALSE;
491     rctx->check          = 0;
492     rctx->oldcapo        = 0;
493     rctx->capo           = 0;
494     rctx->info           = 2;
495     rctx->fine           = (s->stride > 1) ? s->stride : s->total_steps;
496 
497     s->rctx      = rctx;
498     if (s->stype == REVOLVE_ONLINE) rctx->fine = -1;
499 #endif
500   }
501 
502   s->recompute = PETSC_FALSE;
503 
504   ierr = TSGetStages(ts,&numY,PETSC_IGNORE);CHKERRQ(ierr);
505   ierr = StackCreate(PetscObjectComm((PetscObject)ts),s,s->stacksize,numY);CHKERRQ(ierr);
506   PetscFunctionReturn(0);
507 }
508 
509 #undef __FUNCT__
510 #define __FUNCT__ "TSTrajectorySet_Memory"
511 PetscErrorCode TSTrajectorySet_Memory(TSTrajectory tj,TS ts,PetscInt stepnum,PetscReal time,Vec X)
512 {
513   PetscInt       i;
514   Vec            *Y;
515   PetscReal      timeprev;
516   StackElement   e;
517   Stack          *s = (Stack*)tj->data;
518   PetscInt       localstepnum,id;
519   RevolveCTX     *rctx;
520 #ifdef PETSC_HAVE_REVOLVE
521   PetscInt       whattodo,shift;
522 #endif
523   PetscErrorCode ierr;
524 
525   PetscFunctionBegin;
526   if (!s->recompute) { /* use global stepnum in the forward sweep */
527     ierr = TSGetTotalSteps(ts,&stepnum);CHKERRQ(ierr);
528   }
529   /* for consistency */
530   if (!s->recompute && stepnum == 0) ts->ptime_prev = ts->ptime-ts->time_step;
531 
532   if (s->stype == TWO_LEVEL_REVOLVE) { /* two-level mode */
533     localstepnum = stepnum%s->stride;
534     if (stepnum!=s->total_steps && localstepnum==0 && !s->recompute) {
535       id = stepnum/s->stride;
536       if (s->save_stack) {
537         if (stepnum) {
538 #ifdef PETSC_HAVE_REVOLVE
539           PetscPrintf(PETSC_COMM_WORLD,"\x1B[33mDump stack to file\033[0m\n");
540 #endif
541           ierr = StackDumpAll(ts,s,id);CHKERRQ(ierr);
542           s->top = -1; /* reset top */
543 #ifdef PETSC_HAVE_REVOLVE
544           if (s->stype == TWO_LEVEL_REVOLVE) {
545             revolve_reset();
546             revolve_create_offline(s->stride,s->max_cps_ram);
547             rctx = s->rctx;
548             rctx->check = 0;
549             rctx->capo  = 0;
550             rctx->fine  = s->stride;
551           }
552         }
553 #endif
554       } else {
555         ierr = DumpSingle(ts,s,id);CHKERRQ(ierr);
556         PetscPrintf(PETSC_COMM_WORLD,"\x1B[33mDump a single point to file\033[0m\n");
557       }
558     }
559     /* first forward sweep only checkpoints once in each stride */
560     if (!s->recompute && !s->save_stack) PetscFunctionReturn(0);
561   } else if (s->stype == TWO_LEVEL_NOREVOLVE) {
562     localstepnum = stepnum%s->stride;
563     if (stepnum != s->total_steps && stepnum != 0 && localstepnum==0 && !s->recompute) {
564       id = stepnum/s->stride;
565       if (s->save_stack) {
566         ierr = StackDumpAll(ts,s,id);CHKERRQ(ierr);
567         s->top = -1; /* reset top */
568       } else {
569         ierr = DumpSingle(ts,s,id);CHKERRQ(ierr);
570       }
571     }
572   } else {
573     localstepnum = stepnum;
574   }
575 
576   if (s->stype > TWO_LEVEL_NOREVOLVE) {
577 #ifdef PETSC_HAVE_REVOLVE
578     rctx = s->rctx;
579     if (rctx->reverseonestep && stepnum==s->total_steps) { /* intermediate information is ready inside TS, this happens at last time step */
580       //rctx->reverseonestep = PETSC_FALSE;
581       PetscFunctionReturn(0);
582     }
583     if (rctx->stepsleft != 0) { /* advance the solution without checkpointing anything as Revolve requires */
584       rctx->stepsleft--;
585       PetscFunctionReturn(0);
586     }
587 
588     /* let Revolve determine what to do next */
589     shift         = stepnum-localstepnum;
590     rctx->capo    = localstepnum;
591     rctx->oldcapo = rctx->capo;
592     whattodo = revolve_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where);
593     if (s->stype == REVOLVE_ONLINE && whattodo ==7) whattodo = 2;
594     printwhattodo(whattodo,rctx,shift);
595     if (whattodo == -1) SETERRQ(s->comm,PETSC_ERR_LIB,"Error in the Revolve library");
596     if (whattodo == 1) { /* advance some time steps */
597       if (s->stype == REVOLVE_ONLINE && rctx->capo >= s->total_steps-1) {
598         revolve_turn(s->total_steps,&rctx->capo,&rctx->fine);
599         printwhattodo(whattodo,rctx,shift);
600       }
601       rctx->stepsleft = rctx->capo-rctx->oldcapo-1;
602     }
603     if (whattodo == 3 || whattodo == 4) { /* ready for a reverse step */
604       rctx->reverseonestep = PETSC_TRUE;
605     }
606     if (whattodo == 5) { /* restore a checkpoint and ask Revolve what to do next */
607       rctx->oldcapo = rctx->capo;
608       whattodo = revolve_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where); /* must return 1*/
609       printwhattodo(whattodo,rctx,shift);
610       rctx->stepsleft = rctx->capo-rctx->oldcapo;
611     }
612     if (whattodo == 7) { /* save the checkpoint to disk */
613       ierr = DumpSingle(ts,s,rctx->check);CHKERRQ(ierr);
614       rctx->oldcapo = rctx->capo;
615       whattodo = revolve_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where); /* must return 1*/
616       printwhattodo(whattodo,rctx,shift);
617       rctx->stepsleft = rctx->capo-rctx->oldcapo-1;
618     }
619     if (whattodo != 2) {
620       PetscFunctionReturn(0);
621     } else { /* store a checkpoint to RAM and ask Revolve how many time steps to advance next */
622       rctx->oldcapo = rctx->capo;
623       whattodo = revolve_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where); /* must return 1*/
624       printwhattodo(whattodo,rctx,shift);
625       if (s->stype == REVOLVE_ONLINE && rctx->capo >= s->total_steps-1) {
626         revolve_turn(s->total_steps,&rctx->capo,&rctx->fine);
627         printwhattodo(whattodo,rctx,shift);
628       }
629       rctx->stepsleft = rctx->capo-rctx->oldcapo-1;
630     }
631 #endif
632   } else { /* Revolve is not used */
633     if (s->solution_only) {
634       /* skip the last two steps of each stride or the whole interval */
635       if (stepnum >= s->total_steps-1 || s->recompute) PetscFunctionReturn(0); //?
636     } else {
637       /* skip the first and the last steps of each stride or the whole interval */
638       if (localstepnum == 0 || stepnum == s->total_steps) PetscFunctionReturn(0);
639     }
640   }
641 
642   /* checkpoint to memory, rctx->check gives the index in the stack */
643   if (s->stype == REVOLVE_ONLINE && rctx->check != s->top+1) { /* overwrite some non-top checkpoint in the stack*/
644     ierr = StackFind(s,&e,rctx->check);CHKERRQ(ierr);
645     ierr = VecCopy(X,e->X);CHKERRQ(ierr);
646     e->stepnum  = stepnum;
647     e->time     = time;
648     ierr        = TSGetPrevTime(ts,&timeprev);CHKERRQ(ierr);
649     e->timeprev = timeprev;
650   } else { /* push to stack */
651     if (localstepnum < s->top) SETERRQ(s->comm,PETSC_ERR_MEMC,"Illegal modification of a non-top stack element");
652     ierr = PetscCalloc1(1,&e);CHKERRQ(ierr);
653     ierr = VecDuplicate(X,&e->X);CHKERRQ(ierr);
654     ierr = VecCopy(X,e->X);CHKERRQ(ierr);
655     if (s->numY > 0 && !s->solution_only) {
656       ierr = TSGetStages(ts,&s->numY,&Y);CHKERRQ(ierr);
657       ierr = VecDuplicateVecs(Y[0],s->numY,&e->Y);CHKERRQ(ierr);
658       for (i=0;i<s->numY;i++) {
659         ierr = VecCopy(Y[i],e->Y[i]);CHKERRQ(ierr);
660       }
661     }
662     e->stepnum  = stepnum;
663     e->time     = time;
664     /* for consistency */
665     if (stepnum == 0) {
666       e->timeprev = e->time - ts->time_step;
667     } else {
668       ierr = TSGetPrevTime(ts,&timeprev);CHKERRQ(ierr);
669       e->timeprev = timeprev;
670     }
671     ierr = StackPush(s,e);CHKERRQ(ierr);
672   }
673 
674   PetscFunctionReturn(0);
675 }
676 
677 #undef __FUNCT__
678 #define __FUNCT__ "TSTrajectoryGet_Memory"
679 PetscErrorCode TSTrajectoryGet_Memory(TSTrajectory tj,TS ts,PetscInt stepnum,PetscReal *t)
680 {
681   Vec            *Y;
682   PetscInt       i;
683   StackElement   e = NULL;
684   Stack          *s = (Stack*)tj->data;
685   PetscReal      stepsize;
686   PetscInt       localstepnum,id,steps;
687 #ifdef PETSC_HAVE_REVOLVE
688   PetscInt       whattodo,shift;
689 #endif
690   PetscErrorCode ierr;
691 
692   PetscFunctionBegin;
693   ierr = TSGetTotalSteps(ts,&stepnum);CHKERRQ(ierr);
694   if (stepnum == 0) PetscFunctionReturn(0);
695   if (s->stype == TWO_LEVEL_REVOLVE) { /* two-level mode */
696     localstepnum = stepnum%s->stride;
697     if (localstepnum == 0) {
698 #ifdef PETSC_HAVE_REVOLVE
699       PetscPrintf(PETSC_COMM_WORLD,"\x1B[33mLoad stack from file\033[0m\n");
700 #endif
701       if (s->save_stack) {
702         id = stepnum/s->stride;
703         ierr = StackLoadAll(ts,s,id);CHKERRQ(ierr);
704         s->top = s->stacksize-1;
705       } else {
706         /* save ts context */
707         steps = ts->steps;
708         id = stepnum/s->stride-1;
709         ierr = LoadSingle(ts,s,id);CHKERRQ(ierr);
710         ts->steps = ts->total_steps;
711         s->recompute = PETSC_TRUE;
712       }
713 
714 #ifdef PETSC_HAVE_REVOLVE
715       revolve_reset();
716       revolve_create_offline(s->stride,s->max_cps_ram);
717       s->rctx->check = 0;
718       s->rctx->capo  = 0;
719       s->rctx->fine  = s->stride;
720       if (!s->solution_only) {
721         whattodo = 0;
722         while(whattodo!=3) { /* stupid revolve */
723           whattodo = revolve_action(&s->rctx->check,&s->rctx->capo,&s->rctx->fine,s->rctx->snaps_in,&s->rctx->info,&s->rctx->where);
724         }
725       }
726 #endif
727       for (i=ts->steps;i<stepnum;i++) { /* assume fixed step size */
728         ierr = TSTrajectorySet(ts->trajectory,ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
729         ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
730         ierr = TSStep(ts);CHKERRQ(ierr);
731         if (ts->event) {
732           ierr = TSEventMonitor(ts);CHKERRQ(ierr);
733         }
734         if (!ts->steprollback) {
735           ierr = TSPostStep(ts);CHKERRQ(ierr);
736         }
737       }
738       ts->steps = steps;
739       ts->total_steps = stepnum;
740     }
741   } else if (s->stype == TWO_LEVEL_NOREVOLVE) {
742     localstepnum = stepnum%s->stride;
743     if (stepnum != s->total_steps && localstepnum==0) {
744       id = stepnum/s->stride;
745       if (s->save_stack) {
746         ierr = StackLoadAll(ts,s,id);CHKERRQ(ierr);
747       } else {
748         ierr = LoadSingle(ts,s,id);CHKERRQ(ierr);
749       }
750     }
751   } else {
752     localstepnum = stepnum;
753   }
754 #ifdef PETSC_HAVE_REVOLVE
755   if (s->stype > TWO_LEVEL_NOREVOLVE && s->rctx->reverseonestep) { /* ready for the reverse step */
756     ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr);
757     ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr);
758     s->rctx->reverseonestep = PETSC_FALSE;
759     PetscFunctionReturn(0);
760   }
761 #endif
762   if ((!s->solution_only && localstepnum == 0) || stepnum == s->total_steps) {
763     ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr);
764     ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr);
765     PetscFunctionReturn(0);
766   }
767 
768   /* recomputation is needed */
769   steps = ts->steps; /* save the TS context */
770   s->recompute = PETSC_TRUE;
771 
772   if (s->stype > TWO_LEVEL_NOREVOLVE) {
773 #ifdef PETSC_HAVE_REVOLVE
774     s->rctx->capo = localstepnum;
775     shift = stepnum-localstepnum;
776     whattodo = revolve_action(&s->rctx->check,&s->rctx->capo,&s->rctx->fine,s->rctx->snaps_in,&s->rctx->info,&s->rctx->where);
777     if (s->stype == REVOLVE_ONLINE && whattodo == 8) whattodo = 5;
778     printwhattodo(whattodo,s->rctx,shift);
779 #endif
780     if (s->stype == REVOLVE_MULTISTAGE && !s->rctx->where) {
781       ierr = LoadSingle(ts,s,s->rctx->check);CHKERRQ(ierr);
782       ts->steps = ts->total_steps;
783     }
784   }
785   /* restore a checkpoint */
786   if (!(s->stype == REVOLVE_MULTISTAGE && !s->rctx->where)) {
787     if (s->stype == REVOLVE_ONLINE) {
788       ierr = StackFind(s,&e,s->rctx->check);CHKERRQ(ierr);
789     } else {
790       ierr = StackTop(s,&e);CHKERRQ(ierr);
791     }
792     ierr = VecCopy(e->X,ts->vec_sol);CHKERRQ(ierr);
793     if (!s->solution_only) {
794       ierr = TSGetStages(ts,&s->numY,&Y);CHKERRQ(ierr);
795       for (i=0;i<s->numY;i++) {
796         ierr = VecCopy(e->Y[i],Y[i]);CHKERRQ(ierr);
797       }
798     }
799     *t = e->time;
800     ierr = TSSetTimeStep(ts,(*t)-e->timeprev);CHKERRQ(ierr);
801     /* reset ts context */
802     ts->steps      = e->stepnum; /* global stepnum */
803     ts->ptime      = e->time;
804     ts->ptime_prev = e->timeprev;
805   }
806   for (i=ts->steps;i<stepnum;i++) { /* assume fixed step size */
807     ierr = TSTrajectorySet(ts->trajectory,ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
808     ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
809     ierr = TSStep(ts);CHKERRQ(ierr);
810     if (ts->event) {
811       ierr = TSEventMonitor(ts);CHKERRQ(ierr);
812     }
813     if (!ts->steprollback) {
814       ierr = TSPostStep(ts);CHKERRQ(ierr);
815     }
816   }
817   /* reverseonestep must be true after the for loop */
818   ts->steps = steps;
819   ts->total_steps = stepnum;
820   ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr);
821   ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr);
822   if (e && e->stepnum >= stepnum) {
823     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);
824   }
825   if (s->stype != REVOLVE_ONLINE && e && stepnum-e->stepnum==1) { /* the restored checkpoint can be deleted now */
826     ierr = StackPop(s,&e);CHKERRQ(ierr);
827     ierr = VecDestroy(&e->X);CHKERRQ(ierr);
828     if (!s->solution_only) {
829       ierr = VecDestroyVecs(s->numY,&e->Y);CHKERRQ(ierr);
830     }
831     ierr = PetscFree(e);CHKERRQ(ierr);
832   }
833   if (s->stype > TWO_LEVEL_NOREVOLVE) {
834 #ifdef PETSC_HAVE_REVOLVE
835     s->rctx->reverseonestep = PETSC_FALSE;
836 #endif
837   }
838   PetscFunctionReturn(0);
839 }
840 
841 #undef __FUNCT__
842 #define __FUNCT__ "TSTrajectoryDestroy_Memory"
843 PETSC_EXTERN PetscErrorCode TSTrajectoryDestroy_Memory(TSTrajectory tj)
844 {
845   Stack          *s = (Stack*)tj->data;
846   PetscErrorCode ierr;
847 
848   PetscFunctionBegin;
849   if (s->stype > TWO_LEVEL_NOREVOLVE) {
850 #ifdef PETSC_HAVE_REVOLVE
851     revolve_reset();
852 #endif
853   }
854   ierr = StackDestroy(s);CHKERRQ(ierr);
855   PetscFunctionReturn(0);
856 }
857 
858 /*MC
859       TSTRAJECTORYMEMORY - Stores each solution of the ODE/ADE in memory
860 
861   Level: intermediate
862 
863 .seealso:  TSTrajectoryCreate(), TS, TSTrajectorySetType()
864 
865 M*/
866 #undef __FUNCT__
867 #define __FUNCT__ "TSTrajectoryCreate_Memory"
868 PETSC_EXTERN PetscErrorCode TSTrajectoryCreate_Memory(TSTrajectory tj,TS ts)
869 {
870   Stack          *s;
871   PetscErrorCode ierr;
872 
873   PetscFunctionBegin;
874   tj->ops->set            = TSTrajectorySet_Memory;
875   tj->ops->get            = TSTrajectoryGet_Memory;
876   tj->ops->setup          = TSTrajectorySetUp_Memory;
877   tj->ops->destroy        = TSTrajectoryDestroy_Memory;
878   tj->ops->setfromoptions = TSTrajectorySetFromOptions_Memory;
879 
880   ierr = PetscCalloc1(1,&s);CHKERRQ(ierr);
881   s->stype        = NONE;
882   s->max_cps_ram  = -1; /* -1 indicates that it is not set */
883   s->max_cps_disk = -1; /* -1 indicates that it is not set */
884   s->stride       = 0; /* if not zero, two-level checkpointing will be used */
885   s->use_online   = PETSC_FALSE;
886   s->save_stack   = PETSC_TRUE;
887   s->solution_only= PETSC_TRUE;
888   tj->data        = s;
889   PetscFunctionReturn(0);
890 }
891