xref: /petsc/src/ts/trajectory/impls/memory/trajmemory.c (revision 28816579e4dba1d3c21ae89f387694e9bc6e0621)
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     skip_trajectory;
37   PetscBool     solution_only;
38   PetscBool     save_stack;
39   MPI_Comm      comm;
40   RevolveCTX    *rctx;
41   PetscInt      max_cps_ram;  /* maximum checkpoints in RAM */
42   PetscInt      max_cps_disk; /* maximum checkpoints on disk */
43   PetscInt      stride;
44   PetscInt      total_steps;  /* total number of steps */
45   PetscInt      numY;
46   PetscInt      stacksize;
47   PetscInt      top;          /* top of the stack */
48   StackElement  *stack;       /* container */
49 } Stack;
50 
51 #ifdef PETSC_HAVE_REVOLVE
52 static void printwhattodo(PetscInt whattodo,RevolveCTX *rctx,PetscInt shift)
53 {
54   switch(whattodo) {
55     case 1:
56       PetscPrintf(PETSC_COMM_WORLD,"\x1B[35mAdvance from %D to %D\033[0m\n",rctx->oldcapo+shift,rctx->capo+shift);
57       break;
58     case 2:
59       PetscPrintf(PETSC_COMM_WORLD,"\x1B[35mStore in checkpoint number %D (located in RAM)\033[0m\n",rctx->check);
60       break;
61     case 3:
62       PetscPrintf(PETSC_COMM_WORLD,"\x1B[35mFirst turn: Initialize adjoints and reverse first step\033[0m\n");
63       break;
64     case 4:
65       PetscPrintf(PETSC_COMM_WORLD,"\x1B[35mForward and reverse one step\033[0m\n");
66       break;
67     case 5:
68       PetscPrintf(PETSC_COMM_WORLD,"\x1B[35mRestore in checkpoint number %D (located in RAM)\033[0m\n",rctx->check);
69       break;
70     case 7:
71       PetscPrintf(PETSC_COMM_WORLD,"\x1B[35mStore in checkpoint number %D (located on disk)\033[0m\n",rctx->check);
72       break;
73     case 8:
74       PetscPrintf(PETSC_COMM_WORLD,"\x1B[35mRestore in checkpoint number %D (located on disk)\033[0m\n",rctx->check);
75       break;
76     case -1:
77       PetscPrintf(PETSC_COMM_WORLD,"\x1B[35mError!");
78       break;
79   }
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_prev-ts->ptime);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 #ifdef PETSC_HAVE_REVOLVE
435   RevolveCTX     *rctx;
436 #endif
437   PetscInt       numY;
438   PetscBool      flg;
439   PetscErrorCode ierr;
440 
441   PetscFunctionBegin;
442   PetscStrcmp(((PetscObject)ts->adapt)->type_name,TSADAPTNONE,&flg);
443   if (flg) { /* fixed time step */
444     s->total_steps = PetscMin(ts->max_steps,(PetscInt)(ceil(ts->max_time/ts->time_step)));
445   }
446   if (s->max_cps_ram > 1) s->stacksize = s->max_cps_ram;
447   if (s->stride > 1) { /* two level mode works for both fixed time step and adaptive time step */
448     if (s->max_cps_ram > 1 && s->max_cps_ram < s->stride-1) { /* use revolve_offline for each stride */
449       s->stype = TWO_LEVEL_REVOLVE;
450     }else { /* checkpoint all for each stride */
451       s->stype     = TWO_LEVEL_NOREVOLVE;
452       s->stacksize = s->stride-1;
453     }
454   } else {
455     if (flg) { /* fixed time step */
456       if (s->max_cps_ram >= s->total_steps-1 || s->max_cps_ram < 1) { /* checkpoint all */
457         s->stype     = NONE;
458         s->stacksize = s->solution_only ? s->total_steps : s->total_steps-1;
459       } else {
460         if (s->max_cps_disk > 1) { /* disk can be used */
461           s->stype = REVOLVE_MULTISTAGE;
462         } else { /* memory only */
463           s->stype = REVOLVE_OFFLINE;
464         }
465       }
466     } else { /* adaptive time step */
467       s->stype = REVOLVE_ONLINE;
468     }
469     if (s->use_online) { /* trick into online */
470       s->stype     = REVOLVE_ONLINE;
471       s->stacksize = s->max_cps_ram;
472     }
473   }
474 
475   if (s->stype > TWO_LEVEL_NOREVOLVE) {
476 #ifndef PETSC_HAVE_REVOLVE
477     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.");
478 #else
479     if (s->stype == TWO_LEVEL_REVOLVE) revolve_create_offline(s->stride,s->max_cps_ram);
480     else if (s->stype == REVOLVE_OFFLINE) revolve_create_offline(s->total_steps,s->max_cps_ram);
481     else if (s->stype == REVOLVE_ONLINE) revolve_create_online(s->max_cps_ram);
482     else if (s->stype == REVOLVE_MULTISTAGE) revolve_create_multistage(s->total_steps,s->max_cps_ram+s->max_cps_disk,s->max_cps_ram);
483 
484     ierr = PetscCalloc1(1,&rctx);CHKERRQ(ierr);
485     rctx->snaps_in       = s->max_cps_ram; /* for theta methods snaps_in=2*max_cps_ram */
486     rctx->reverseonestep = PETSC_FALSE;
487     rctx->check          = 0;
488     rctx->oldcapo        = 0;
489     rctx->capo           = 0;
490     rctx->info           = 2;
491     rctx->fine           = (s->stride > 1) ? s->stride : s->total_steps;
492 
493     s->rctx      = rctx;
494     if (s->stype == REVOLVE_ONLINE) rctx->fine = -1;
495 #endif
496   }
497 
498   s->recompute = PETSC_FALSE;
499 
500   ierr = TSGetStages(ts,&numY,PETSC_IGNORE);CHKERRQ(ierr);
501   ierr = StackCreate(PetscObjectComm((PetscObject)ts),s,s->stacksize,numY);CHKERRQ(ierr);
502   PetscFunctionReturn(0);
503 }
504 
505 #undef __FUNCT__
506 #define __FUNCT__ "ApplyRevolve"
507 static PetscErrorCode ApplyRevolve(Stack *s,PetscInt stepnum,PetscInt localstepnum,PetscInt *store)
508 {
509 #ifdef PETSC_HAVE_REVOLVE
510   PetscInt       shift,whattodo;
511   RevolveCTX     *rctx;
512 #endif
513 
514   PetscFunctionBegin;
515   *store = 0;
516 #ifdef PETSC_HAVE_REVOLVE
517   rctx = s->rctx;
518   if (rctx->reverseonestep && stepnum == s->total_steps) { /* intermediate information is ready inside TS, this happens at last time step */
519     rctx->reverseonestep = PETSC_FALSE;
520     PetscFunctionReturn(0);
521   }
522   if (rctx->stepsleft > 0) { /* advance the solution without checkpointing anything as Revolve requires */
523     rctx->stepsleft--;
524     PetscFunctionReturn(0);
525   }
526   /* let Revolve determine what to do next */
527   shift         = stepnum-localstepnum;
528   rctx->oldcapo = rctx->capo;
529   rctx->capo    = localstepnum;
530   whattodo = revolve_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where);
531   if (s->stype == REVOLVE_ONLINE && whattodo == 8) whattodo = 5;
532   if (s->stype == REVOLVE_ONLINE && whattodo == 7) whattodo = 2;
533   printwhattodo(whattodo,rctx,shift);
534   if (whattodo == -1) SETERRQ(s->comm,PETSC_ERR_LIB,"Error in the Revolve library");
535   if (whattodo == 1) { /* advance some time steps */
536     if (s->stype == REVOLVE_ONLINE && rctx->capo >= s->total_steps-1) {
537       revolve_turn(s->total_steps,&rctx->capo,&rctx->fine);
538       printwhattodo(whattodo,rctx,shift);
539     }
540     rctx->stepsleft = rctx->capo-rctx->oldcapo-1;
541   }
542   if (whattodo == 3 || whattodo == 4) { /* ready for a reverse step */
543     rctx->reverseonestep = PETSC_TRUE;
544   }
545   if (whattodo == 5) { /* restore a checkpoint and ask Revolve what to do next */
546     rctx->oldcapo = rctx->capo;
547     whattodo = revolve_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where); /* must return 1 or 3 or 4*/
548     printwhattodo(whattodo,rctx,shift);
549     if (whattodo == 3 || whattodo == 4) rctx->reverseonestep = PETSC_TRUE;
550     if (whattodo == 1) rctx->stepsleft = rctx->capo-rctx->oldcapo;
551   }
552   if (whattodo == 7) { /* save the checkpoint to disk */
553     *store = 2;
554     rctx->oldcapo = rctx->capo;
555     whattodo = revolve_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where); /* must return 1*/
556     printwhattodo(whattodo,rctx,shift);
557     rctx->stepsleft = rctx->capo-rctx->oldcapo-1;
558   }
559   if (whattodo == 2) { /* store a checkpoint to RAM and ask Revolve how many time steps to advance next */
560     *store = 1;
561     rctx->oldcapo = rctx->capo;
562     whattodo = revolve_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where); /* must return 1*/
563     printwhattodo(whattodo,rctx,shift);
564     if (s->stype == REVOLVE_ONLINE && rctx->capo >= s->total_steps-1) {
565       revolve_turn(s->total_steps,&rctx->capo,&rctx->fine);
566       printwhattodo(whattodo,rctx,shift);
567     }
568     rctx->stepsleft = rctx->capo-rctx->oldcapo-1;
569   }
570 #endif
571   PetscFunctionReturn(0);
572 }
573 
574 #undef __FUNCT__
575 #define __FUNCT__ "ElementCreate"
576 static PetscErrorCode ElementCreate(TS ts,Stack *s,StackElement *e,PetscInt stepnum,PetscReal time,Vec X)
577 {
578   Vec            *Y;
579   PetscInt       i;
580   PetscReal      timeprev;
581   PetscErrorCode ierr;
582 
583   PetscFunctionBegin;
584   ierr = PetscCalloc1(1,e);CHKERRQ(ierr);
585   ierr = VecDuplicate(X,&(*e)->X);CHKERRQ(ierr);
586   ierr = VecCopy(X,(*e)->X);CHKERRQ(ierr);
587   if (s->numY > 0 && !s->solution_only) {
588     ierr = TSGetStages(ts,&s->numY,&Y);CHKERRQ(ierr);
589     ierr = VecDuplicateVecs(Y[0],s->numY,&(*e)->Y);CHKERRQ(ierr);
590     for (i=0;i<s->numY;i++) {
591       ierr = VecCopy(Y[i],(*e)->Y[i]);CHKERRQ(ierr);
592     }
593   }
594   (*e)->stepnum = stepnum;
595   (*e)->time    = time;
596   /* for consistency */
597   if (stepnum == 0) {
598     (*e)->timeprev = (*e)->time - ts->time_step;
599   } else {
600     ierr = TSGetPrevTime(ts,&timeprev);CHKERRQ(ierr);
601     (*e)->timeprev = timeprev;
602   }
603   PetscFunctionReturn(0);
604 }
605 
606 #undef __FUNCT__
607 #define __FUNCT__ "ElementDestroy"
608 static PetscErrorCode ElementDestroy(Stack *s,StackElement e)
609 {
610   PetscErrorCode ierr;
611 
612   PetscFunctionBegin;
613   ierr = VecDestroy(&e->X);CHKERRQ(ierr);
614   if (!s->solution_only) {
615     ierr = VecDestroyVecs(s->numY,&e->Y);CHKERRQ(ierr);
616   }
617   ierr = PetscFree(e);CHKERRQ(ierr);
618   PetscFunctionReturn(0);
619 }
620 
621 #undef __FUNCT__
622 #define __FUNCT__ "UpdateTS"
623 static PetscErrorCode UpdateTS(TS ts,Stack *s,StackElement e)
624 {
625   Vec            *Y;
626   PetscInt       i;
627   PetscErrorCode ierr;
628 
629   PetscFunctionBegin;
630   ierr = VecCopy(e->X,ts->vec_sol);CHKERRQ(ierr);
631   if (!s->solution_only) {
632     ierr = TSGetStages(ts,&s->numY,&Y);CHKERRQ(ierr);
633     for (i=0;i<s->numY;i++) {
634       ierr = VecCopy(e->Y[i],Y[i]);CHKERRQ(ierr);
635     }
636   }
637   ierr = TSSetTimeStep(ts,e->timeprev-e->time);CHKERRQ(ierr); /* stepsize will be negative */
638   ts->ptime      = e->time;
639   ts->ptime_prev = e->timeprev;
640   PetscFunctionReturn(0);
641 }
642 
643 #undef __FUNCT__
644 #define __FUNCT__ "ReCompute"
645 static PetscErrorCode ReCompute(TS ts,Stack *s,PetscInt stepnumbegin,PetscInt stepnumend)
646 {
647   PetscInt       i,adjsteps;
648   PetscReal      stepsize;
649   PetscErrorCode ierr;
650 
651   PetscFunctionBegin;
652   adjsteps = ts->steps;
653   /* reset ts context */
654   ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr);
655   ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr);
656   ts->steps = stepnumbegin; /* global step number */
657   for (i=ts->steps;i<stepnumend;i++) { /* assume fixed step size */
658     if (s->solution_only && !s->skip_trajectory) { /* revolve online need this */
659       ierr = TSTrajectorySet(ts->trajectory,ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
660     }
661     ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
662     ierr = TSStep(ts);CHKERRQ(ierr);
663     if (!s->solution_only && !s->skip_trajectory) {
664       ierr = TSTrajectorySet(ts->trajectory,ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
665     }
666     if (ts->event) {
667       ierr = TSEventMonitor(ts);CHKERRQ(ierr);
668     }
669     if (!ts->steprollback) {
670       ierr = TSPostStep(ts);CHKERRQ(ierr);
671     }
672   }
673   ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr);
674   ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr);
675   ts->steps = adjsteps;
676   ts->total_steps = stepnumend;
677   PetscFunctionReturn(0);
678 }
679 
680 #undef __FUNCT__
681 #define __FUNCT__ "SetTrajN"
682 static PetscErrorCode SetTrajN(TS ts,Stack *s,PetscInt stepnum,PetscReal time,Vec X)
683 {
684   StackElement   e;
685   PetscErrorCode ierr;
686 
687   PetscFunctionBegin;
688   if (s->solution_only) {
689     /* skip the last two steps of each stride or the whole interval */
690     if (stepnum >= s->total_steps-1 || s->recompute) PetscFunctionReturn(0); //?
691   } else {
692     /* skip the first and the last steps of each stride or the whole interval */
693     if (stepnum == 0 || stepnum == s->total_steps) PetscFunctionReturn(0);
694   }
695   if (stepnum < s->top) SETERRQ(s->comm,PETSC_ERR_MEMC,"Illegal modification of a non-top stack element");
696   ierr = ElementCreate(ts,s,&e,stepnum,time,X);CHKERRQ(ierr);
697   ierr = StackPush(s,e);CHKERRQ(ierr);
698   PetscFunctionReturn(0);
699 }
700 
701 #undef __FUNCT__
702 #define __FUNCT__ "GetTrajN"
703 static PetscErrorCode GetTrajN(TS ts,Stack *s,PetscInt stepnum)
704 {
705   PetscReal      stepsize;
706   StackElement   e;
707   PetscErrorCode ierr;
708 
709   PetscFunctionBegin;
710   if (stepnum == s->total_steps) {
711     ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr);
712     ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr);
713     PetscFunctionReturn(0);
714   }
715   /* restore a checkpoint */
716   ierr = StackTop(s,&e);CHKERRQ(ierr);
717   ierr = UpdateTS(ts,s,e);CHKERRQ(ierr);
718   if (s->solution_only) {/* recompute one step */
719     s->recompute = PETSC_TRUE;
720     ierr = ReCompute(ts,s,e->stepnum,stepnum);CHKERRQ(ierr);
721   }
722   ierr = StackPop(s,&e);CHKERRQ(ierr);
723   ierr = ElementDestroy(s,e);CHKERRQ(ierr);
724   PetscFunctionReturn(0);
725 }
726 
727 #undef __FUNCT__
728 #define __FUNCT__ "SetTrajROF"
729 static PetscErrorCode SetTrajROF(TS ts,Stack *s,PetscInt stepnum,PetscReal time,Vec X)
730 {
731   PetscInt       store;
732   StackElement   e;
733   PetscErrorCode ierr;
734 
735   PetscFunctionBegin;
736   if (!s->solution_only && stepnum == 0) PetscFunctionReturn(0);
737   ierr = ApplyRevolve(s,stepnum,stepnum,&store);CHKERRQ(ierr);
738   if (store == 1) {
739     if (stepnum < s->top) SETERRQ(s->comm,PETSC_ERR_MEMC,"Illegal modification of a non-top stack element");
740     ierr = ElementCreate(ts,s,&e,stepnum,time,X);CHKERRQ(ierr);
741     ierr = StackPush(s,e);CHKERRQ(ierr);
742   }
743   PetscFunctionReturn(0);
744 }
745 
746 #undef __FUNCT__
747 #define __FUNCT__ "GetTrajROF"
748 static PetscErrorCode GetTrajROF(TS ts,Stack *s,PetscInt stepnum)
749 {
750 #ifdef PETSC_HAVE_REVOLVE
751   PetscInt       whattodo,shift,store;
752 #endif
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,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,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,laststridesize;
895   StackElement   e;
896 #ifdef PETSC_HAVE_REVOLVE
897   RevolveCTX     *rctx = s->rctx;
898   PetscBool      resetrevolve = PETSC_FALSE;
899 #endif
900 
901   PetscErrorCode ierr;
902 
903   PetscFunctionBegin;
904   localstepnum = stepnum%s->stride;
905   id           = stepnum/s->stride; /* stride index */
906 
907   /* (stride size-1) checkpoints are saved in each stride */
908   laststridesize = s->total_steps%s->stride;
909   if (!laststridesize) laststridesize = s->stride;
910   if (s->solution_only) {
911     if (stepnum == s->total_steps) PetscFunctionReturn(0);
912     if (s->save_stack) {
913       if (!s->recompute && localstepnum == s->stride-1 && stepnum < s->total_steps-laststridesize) {
914 #ifdef PETSC_HAVE_REVOLVE
915         PetscPrintf(PETSC_COMM_WORLD,"\x1B[33mDump stack to file\033[0m\n");
916         resetrevolve = PETSC_TRUE;
917 #endif
918         ierr = StackDumpAll(ts,s,id+1);CHKERRQ(ierr);
919       }
920     } else {
921       if (!s->recompute && localstepnum == 0 && stepnum < s->total_steps-laststridesize ) {
922         PetscPrintf(PETSC_COMM_WORLD,"\x1B[33mDump a single point to file\033[0m\n");
923         ierr = DumpSingle(ts,s,id+1);CHKERRQ(ierr);
924       }
925       if (stepnum < s->total_steps-laststridesize && !s->recompute) PetscFunctionReturn(0); /* no need to checkpoint except last stride in the first sweep */
926     }
927   } else {
928     if (stepnum == 0) PetscFunctionReturn(0);
929     if (s->save_stack) {
930       if (!s->recompute && localstepnum == 0 && stepnum != s->total_steps) { /* do not dump stack for last stride */
931 #ifdef PETSC_HAVE_REVOLVE
932         PetscPrintf(PETSC_COMM_WORLD,"\x1B[33mDump stack to file\033[0m\n");
933         resetrevolve = PETSC_TRUE;
934 #endif
935         ierr = StackDumpAll(ts,s,id);CHKERRQ(ierr);
936       }
937     } else {
938       if (!s->recompute && localstepnum == 1 && stepnum <  s->total_steps-laststridesize ) { /* skip last stride */
939         PetscPrintf(PETSC_COMM_WORLD,"\x1B[33mDump a single point to file\033[0m\n");
940         ierr = DumpSingle(ts,s,id);CHKERRQ(ierr);
941       }
942       if (stepnum <= s->total_steps-laststridesize && !s->recompute) PetscFunctionReturn(0);
943     }
944   }
945 
946   ierr = ApplyRevolve(s,stepnum,localstepnum,&store);CHKERRQ(ierr);
947   if (store == 1){
948     if (localstepnum < s->top) SETERRQ(s->comm,PETSC_ERR_MEMC,"Illegal modification of a non-top stack element");
949     ierr = ElementCreate(ts,s,&e,stepnum,time,X);CHKERRQ(ierr);
950     ierr = StackPush(s,e);CHKERRQ(ierr);
951   }
952 #ifdef PETSC_HAVE_REVOLVE
953   if (resetrevolve) {
954     revolve_reset();
955     revolve_create_offline(s->stride,s->max_cps_ram);
956     rctx = s->rctx;
957     rctx->check = 0;
958     rctx->capo  = 0;
959     rctx->fine  = s->stride;
960     if (s->rctx->reverseonestep) s->rctx->reverseonestep = PETSC_FALSE;
961   }
962 #endif
963   PetscFunctionReturn(0);
964 }
965 
966 #undef __FUNCT__
967 #define __FUNCT__ "GetTrajTLR"
968 static PetscErrorCode GetTrajTLR(TS ts,Stack *s,PetscInt stepnum)
969 {
970 #ifdef PETSC_HAVE_REVOLVE
971   PetscInt       whattodo,shift;
972 #endif
973   PetscInt       localstepnum,id,laststridesize,store;
974   PetscReal      stepsize;
975   StackElement   e;
976   PetscErrorCode ierr;
977 
978   PetscFunctionBegin;
979   localstepnum = stepnum%s->stride;
980   id           = stepnum/s->stride;
981   if (stepnum == s->total_steps) {
982     ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr);
983     ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr);
984 #ifdef PETSC_HAVE_REVOLVE
985     if ( s->rctx->reverseonestep) s->rctx->reverseonestep = PETSC_FALSE;
986 #endif
987     PetscFunctionReturn(0);
988   }
989   laststridesize = s->total_steps%s->stride;
990   if (!laststridesize) laststridesize = s->stride;
991   if (s->solution_only) {
992     /* fill stack */
993     if (localstepnum == 0 && stepnum <= s->total_steps-laststridesize) {
994       if (s->save_stack) {
995 #ifdef PETSC_HAVE_REVOLVE
996         PetscPrintf(PETSC_COMM_WORLD,"\x1B[33mLoad stack from file\033[0m\n");
997 #endif
998         ierr = StackLoadAll(ts,s,id);CHKERRQ(ierr);
999 #ifdef PETSC_HAVE_REVOLVE
1000         revolve_reset();
1001         revolve_create_offline(s->stride,s->max_cps_ram);
1002         s->rctx->check = 0;
1003         s->rctx->capo  = 0;
1004         s->rctx->fine  = s->stride;
1005         whattodo = 0;
1006         while(whattodo!=3) { /* stupid revolve */
1007           whattodo = revolve_action(&s->rctx->check,&s->rctx->capo,&s->rctx->fine,s->rctx->snaps_in,&s->rctx->info,&s->rctx->where);
1008         }
1009 #endif
1010         s->recompute = PETSC_TRUE;
1011         s->skip_trajectory = PETSC_TRUE;
1012         ierr = ReCompute(ts,s,id*s->stride-1,id*s->stride);CHKERRQ(ierr);
1013         s->skip_trajectory = PETSC_FALSE;
1014       } else {
1015 #ifdef PETSC_HAVE_REVOLVE
1016         PetscPrintf(PETSC_COMM_WORLD,"\x1B[33mLoad a single point from file\033[0m\n");
1017 #endif
1018         ierr = LoadSingle(ts,s,id);CHKERRQ(ierr);
1019 #ifdef PETSC_HAVE_REVOLVE
1020         revolve_reset();
1021         revolve_create_offline(s->stride,s->max_cps_ram);
1022         s->rctx->check = 0;
1023         s->rctx->capo  = 0;
1024         s->rctx->fine  = s->stride;
1025 #endif
1026         s->recompute = PETSC_TRUE;
1027         ierr = ReCompute(ts,s,(id-1)*s->stride,id*s->stride);CHKERRQ(ierr);
1028       }
1029       PetscFunctionReturn(0);
1030     }
1031     /* restore a checkpoint */
1032     ierr = StackTop(s,&e);CHKERRQ(ierr);
1033     ierr = UpdateTS(ts,s,e);CHKERRQ(ierr);
1034 #ifdef PETSC_HAVE_REVOLVE
1035     /* start with restoring a checkpoint */
1036     s->rctx->capo = stepnum;
1037     s->rctx->oldcapo = s->rctx->capo;
1038     shift = stepnum-localstepnum;
1039     whattodo = revolve_action(&s->rctx->check,&s->rctx->capo,&s->rctx->fine,s->rctx->snaps_in,&s->rctx->info,&s->rctx->where);
1040     printwhattodo(whattodo,s->rctx,shift);
1041 #endif
1042     s->recompute = PETSC_TRUE;
1043     ierr = ReCompute(ts,s,e->stepnum,stepnum);CHKERRQ(ierr);
1044     if (e->stepnum+1 == stepnum) {
1045       ierr = StackPop(s,&e);CHKERRQ(ierr);
1046       ierr = ElementDestroy(s,e);CHKERRQ(ierr);
1047     }
1048   } else {
1049     /* fill stack with info */
1050     if (localstepnum == 0 && s->total_steps-stepnum >= laststridesize) {
1051       if (s->save_stack) {
1052 #ifdef PETSC_HAVE_REVOLVE
1053         PetscPrintf(PETSC_COMM_WORLD,"\x1B[33mLoad stack from file\033[0m\n");
1054 #endif
1055         ierr = StackLoadAll(ts,s,id);CHKERRQ(ierr);
1056 #ifdef PETSC_HAVE_REVOLVE
1057         revolve_reset();
1058         revolve_create_offline(s->stride,s->max_cps_ram);
1059         s->rctx->check = 0;
1060         s->rctx->capo  = 0;
1061         s->rctx->fine  = s->stride;
1062         whattodo = 0;
1063         while(whattodo!=3) { /* stupid revolve */
1064           whattodo = revolve_action(&s->rctx->check,&s->rctx->capo,&s->rctx->fine,s->rctx->snaps_in,&s->rctx->info,&s->rctx->where);
1065         }
1066 #endif
1067       } else {
1068 #ifdef PETSC_HAVE_REVOLVE
1069         PetscPrintf(PETSC_COMM_WORLD,"\x1B[33mLoad a single point from file\033[0m\n");
1070 #endif
1071         ierr = LoadSingle(ts,s,id-1);CHKERRQ(ierr);
1072 #ifdef PETSC_HAVE_REVOLVE
1073         revolve_reset();
1074         revolve_create_offline(s->stride,s->max_cps_ram);
1075         s->rctx->check = 0;
1076         s->rctx->capo  = 0;
1077         s->rctx->fine  = s->stride;
1078         shift = stepnum-localstepnum;
1079         whattodo = revolve_action(&s->rctx->check,&s->rctx->capo,&s->rctx->fine,s->rctx->snaps_in,&s->rctx->info,&s->rctx->where);
1080         printwhattodo(whattodo,s->rctx,shift);
1081 #endif
1082         ierr = ElementCreate(ts,s,&e,(id-1)*s->stride+1,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
1083         ierr = StackPush(s,e);CHKERRQ(ierr);
1084         s->recompute = PETSC_TRUE;
1085         ierr = ReCompute(ts,s,e->stepnum,id*s->stride);CHKERRQ(ierr);
1086 #ifdef PETSC_HAVE_REVOLVE
1087         if ( s->rctx->reverseonestep) { /* ready for the reverse step */
1088           ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr);
1089           ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr);
1090           s->rctx->reverseonestep = PETSC_FALSE;
1091         }
1092 #endif
1093       }
1094       PetscFunctionReturn(0);
1095     }
1096     /* restore a checkpoint */
1097     ierr = StackTop(s,&e);CHKERRQ(ierr);
1098     ierr = UpdateTS(ts,s,e);CHKERRQ(ierr);
1099     /* 2 revolve actions: restore a checkpoint and then advance */
1100     ierr = ApplyRevolve(s,stepnum,localstepnum,&store);CHKERRQ(ierr);
1101     if (!s->rctx->reverseonestep && s->rctx->stepsleft > 0) {
1102       s->rctx->stepsleft--;
1103       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);
1104     }
1105     if (e->stepnum < stepnum) {
1106       s->recompute = PETSC_TRUE;
1107       ierr = ReCompute(ts,s,e->stepnum,stepnum);CHKERRQ(ierr);
1108     }
1109     if (e->stepnum == stepnum) {
1110       ierr = StackPop(s,&e);CHKERRQ(ierr);
1111       ierr = ElementDestroy(s,e);CHKERRQ(ierr);
1112     }
1113   }
1114 #ifdef PETSC_HAVE_REVOLVE
1115   if (s->rctx->reverseonestep) s->rctx->reverseonestep = PETSC_FALSE;
1116 #endif
1117   PetscFunctionReturn(0);
1118 }
1119 
1120 #undef __FUNCT__
1121 #define __FUNCT__ "SetTrajTLNR"
1122 static PetscErrorCode SetTrajTLNR(TS ts,Stack *s,PetscInt stepnum,PetscReal time,Vec X)
1123 {
1124   PetscInt       localstepnum,id,laststridesize;
1125   StackElement   e;
1126   PetscErrorCode ierr;
1127 
1128   PetscFunctionBegin;
1129   localstepnum = stepnum%s->stride;
1130   id           = stepnum/s->stride;
1131   if (stepnum == s->total_steps) PetscFunctionReturn(0);
1132 
1133   /* (stride size-1) checkpoints are saved in each stride */
1134   laststridesize = s->total_steps%s->stride;
1135   if (!laststridesize) laststridesize = s->stride;
1136   if (s->solution_only) {
1137     if (s->save_stack) {
1138       if (s->recompute) PetscFunctionReturn(0);
1139       if (localstepnum == s->stride-1 && stepnum < s->total_steps-laststridesize) { /* dump when stack is full */
1140         ierr = StackDumpAll(ts,s,id+1);CHKERRQ(ierr);
1141         PetscFunctionReturn(0);
1142       }
1143       if (stepnum == s->total_steps-1) PetscFunctionReturn(0); /* do not checkpoint s->total_steps-1 */
1144     } else {
1145       if (localstepnum == s->stride-1) PetscFunctionReturn(0);
1146       if (!s->recompute && localstepnum == 0 && stepnum < s->total_steps-laststridesize ) {
1147         ierr = DumpSingle(ts,s,id+1);CHKERRQ(ierr);
1148       }
1149       if (stepnum < s->total_steps-laststridesize && !s->recompute) PetscFunctionReturn(0);
1150     }
1151   } else {
1152     if (stepnum == 0) PetscFunctionReturn(0);
1153     if (s->save_stack) {
1154       if (s->recompute) PetscFunctionReturn(0);
1155       if (localstepnum == 0 && stepnum != 0) { /* no stack at point 0 */
1156         ierr = StackDumpAll(ts,s,id);CHKERRQ(ierr);
1157         PetscFunctionReturn(0);
1158       }
1159     } else {
1160       if (localstepnum == 0) PetscFunctionReturn(0); /* skip last point in each stride */
1161       if (!s->recompute && localstepnum == 1 && stepnum < s->total_steps-laststridesize ) { /* skip last stride */
1162         ierr = DumpSingle(ts,s,id);CHKERRQ(ierr);
1163       }
1164       if (stepnum <= s->total_steps-laststridesize && !s->recompute) PetscFunctionReturn(0);
1165     }
1166   }
1167   ierr = ElementCreate(ts,s,&e,stepnum,time,X);CHKERRQ(ierr);
1168   ierr = StackPush(s,e);CHKERRQ(ierr);
1169   PetscFunctionReturn(0);
1170 }
1171 
1172 #undef __FUNCT__
1173 #define __FUNCT__ "GetTrajTLNR"
1174 static PetscErrorCode GetTrajTLNR(TS ts,Stack *s,PetscInt stepnum)
1175 {
1176   PetscInt       id,localstepnum,laststridesize;
1177   PetscReal      stepsize;
1178   StackElement   e;
1179   PetscErrorCode ierr;
1180 
1181   PetscFunctionBegin;
1182   localstepnum = stepnum%s->stride;
1183   if (stepnum == s->total_steps) {
1184     ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr);
1185     ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr);
1186     PetscFunctionReturn(0);
1187   }
1188   laststridesize = s->total_steps%s->stride;
1189   if (!laststridesize) laststridesize = s->stride;
1190   if (s->solution_only) {
1191     /* fill stack with info */
1192     if (localstepnum == 0 && s->total_steps-stepnum >= laststridesize) {
1193       id = stepnum/s->stride;
1194       if (s->save_stack) {
1195         ierr = StackLoadAll(ts,s,id);CHKERRQ(ierr);
1196         s->recompute = PETSC_TRUE;
1197         s->skip_trajectory = PETSC_TRUE;
1198         ierr = ReCompute(ts,s,id*s->stride-1,id*s->stride);CHKERRQ(ierr);
1199         s->skip_trajectory = PETSC_FALSE;
1200       } else {
1201         ierr = LoadSingle(ts,s,id);CHKERRQ(ierr);
1202         s->recompute = PETSC_TRUE;
1203         ierr = ReCompute(ts,s,(id-1)*s->stride,id*s->stride);CHKERRQ(ierr);
1204       }
1205       PetscFunctionReturn(0);
1206     }
1207     /* restore a checkpoint */
1208     ierr = StackPop(s,&e);CHKERRQ(ierr);
1209     ierr = UpdateTS(ts,s,e);CHKERRQ(ierr);
1210     s->recompute = PETSC_TRUE;
1211     s->skip_trajectory = PETSC_TRUE;
1212     ierr = ReCompute(ts,s,e->stepnum,stepnum);CHKERRQ(ierr);
1213     s->skip_trajectory = PETSC_FALSE;
1214     ierr = ElementDestroy(s,e);CHKERRQ(ierr);
1215   } else {
1216     /* fill stack with info */
1217     if (localstepnum == 0 && s->total_steps-stepnum >= laststridesize) {
1218       id = stepnum/s->stride;
1219       if (s->save_stack) {
1220         ierr = StackLoadAll(ts,s,id);CHKERRQ(ierr);
1221       } else {
1222         ierr = LoadSingle(ts,s,id-1);CHKERRQ(ierr);
1223         ierr = ElementCreate(ts,s,&e,(id-1)*s->stride+1,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
1224         ierr = StackPush(s,e);CHKERRQ(ierr);
1225         s->recompute = PETSC_TRUE;
1226         ierr = ReCompute(ts,s,e->stepnum,id*s->stride);CHKERRQ(ierr);
1227       }
1228       PetscFunctionReturn(0);
1229     }
1230     /* restore a checkpoint */
1231     ierr = StackPop(s,&e);CHKERRQ(ierr);
1232     ierr = UpdateTS(ts,s,e);CHKERRQ(ierr);
1233     ierr = ElementDestroy(s,e);CHKERRQ(ierr);
1234   }
1235   PetscFunctionReturn(0);
1236 }
1237 
1238 #undef __FUNCT__
1239 #define __FUNCT__ "SetTrajRMS"
1240 static PetscErrorCode SetTrajRMS(TS ts,Stack *s,PetscInt stepnum,PetscReal time,Vec X)
1241 {
1242   PetscInt       store;
1243   StackElement   e;
1244   PetscErrorCode ierr;
1245 
1246   PetscFunctionBegin;
1247   if (!s->solution_only && stepnum == 0) PetscFunctionReturn(0);
1248   ierr = ApplyRevolve(s,stepnum,stepnum,&store);CHKERRQ(ierr);
1249   if (store == 1){
1250     if (stepnum < s->top) SETERRQ(s->comm,PETSC_ERR_MEMC,"Illegal modification of a non-top stack element");
1251     ierr = ElementCreate(ts,s,&e,stepnum,time,X);CHKERRQ(ierr);
1252     ierr = StackPush(s,e);CHKERRQ(ierr);
1253   } else if (store == 2) {
1254     ierr = DumpSingle(ts,s,s->rctx->check);CHKERRQ(ierr);
1255   }
1256   PetscFunctionReturn(0);
1257 }
1258 
1259 #undef __FUNCT__
1260 #define __FUNCT__ "GetTrajRMS"
1261 static PetscErrorCode GetTrajRMS(TS ts,Stack *s,PetscInt stepnum)
1262 {
1263 #ifdef PETSC_HAVE_REVOLVE
1264   PetscInt       whattodo,shift;
1265 #endif
1266   PetscInt       restart;
1267   PetscBool      ondisk;
1268   PetscReal      stepsize;
1269   StackElement   e;
1270   PetscErrorCode ierr;
1271 
1272   PetscFunctionBegin;
1273   if (stepnum == 0 || stepnum == s->total_steps) {
1274     ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr);
1275     ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr);
1276 #ifdef PETSC_HAVE_REVOLVE
1277     if (s->rctx->reverseonestep) s->rctx->reverseonestep = PETSC_FALSE;
1278 #endif
1279     PetscFunctionReturn(0);
1280   }
1281 #ifdef PETSC_HAVE_REVOLVE
1282   s->rctx->capo = stepnum;
1283   s->rctx->oldcapo = s->rctx->capo;
1284   shift = 0;
1285   whattodo = revolve_action(&s->rctx->check,&s->rctx->capo,&s->rctx->fine,s->rctx->snaps_in,&s->rctx->info,&s->rctx->where); /* whattodo=restore */
1286   printwhattodo(whattodo,s->rctx,shift);
1287 #endif
1288   /* restore a checkpoint */
1289   restart = s->rctx->capo;
1290   if (!s->rctx->where) {
1291     ondisk = PETSC_TRUE;
1292     ierr = LoadSingle(ts,s,s->rctx->check);CHKERRQ(ierr);
1293     ierr = TSSetTimeStep(ts,ts->ptime_prev-ts->ptime);CHKERRQ(ierr);
1294   } else {
1295     ondisk = PETSC_FALSE;
1296     ierr = StackTop(s,&e);CHKERRQ(ierr);
1297     ierr = UpdateTS(ts,s,e);CHKERRQ(ierr);
1298   }
1299 #ifdef PETSC_HAVE_REVOLVE
1300   if (!s->solution_only) { /* whattodo must be 5 or 8 */
1301     /* ask Revolve what to do next */
1302     s->rctx->oldcapo = s->rctx->capo;
1303     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*/
1304     printwhattodo(whattodo,s->rctx,shift);
1305     if (whattodo == 3 || whattodo == 4) s->rctx->reverseonestep = PETSC_TRUE;
1306     if (whattodo == 1) s->rctx->stepsleft = s->rctx->capo-s->rctx->oldcapo;
1307     if (!s->rctx->reverseonestep && s->rctx->stepsleft > 0) {
1308       s->rctx->stepsleft--;
1309       PetscPrintf(PETSC_COMM_WORLD,"\x1B[35mSkip the step from %D to %D (already checkpointed)\033[0m\n",s->rctx->oldcapo,s->rctx->oldcapo+1);
1310     }
1311     restart++; /* skip one step */
1312   }
1313 #endif
1314   if (s->solution_only || (!s->solution_only && restart < stepnum)) {
1315     s->recompute = PETSC_TRUE;
1316     ierr = ReCompute(ts,s,restart,stepnum);CHKERRQ(ierr);
1317   }
1318   if (!ondisk && ( (s->solution_only && e->stepnum+1 == stepnum) || (!s->solution_only && e->stepnum == stepnum) )) {
1319     ierr = StackPop(s,&e);CHKERRQ(ierr);
1320     ierr = ElementDestroy(s,e);CHKERRQ(ierr);
1321   }
1322 #ifdef PETSC_HAVE_REVOLVE
1323   if (s->rctx->reverseonestep) s->rctx->reverseonestep = PETSC_FALSE;
1324 #endif
1325   PetscFunctionReturn(0);
1326 }
1327 
1328 #undef __FUNCT__
1329 #define __FUNCT__ "TSTrajectorySet_Memory"
1330 PetscErrorCode TSTrajectorySet_Memory(TSTrajectory tj,TS ts,PetscInt stepnum,PetscReal time,Vec X)
1331 {
1332   Stack          *s = (Stack*)tj->data;
1333   PetscErrorCode ierr;
1334 
1335   PetscFunctionBegin;
1336   if (!s->recompute) { /* use global stepnum in the forward sweep */
1337     ierr = TSGetTotalSteps(ts,&stepnum);CHKERRQ(ierr);
1338   }
1339   /* for consistency */
1340   if (!s->recompute && stepnum == 0) ts->ptime_prev = ts->ptime-ts->time_step;
1341   switch (s->stype) {
1342     case NONE:
1343       ierr = SetTrajN(ts,s,stepnum,time,X);CHKERRQ(ierr);
1344       break;
1345     case TWO_LEVEL_NOREVOLVE:
1346       ierr = SetTrajTLNR(ts,s,stepnum,time,X);CHKERRQ(ierr);
1347       break;
1348     case TWO_LEVEL_REVOLVE:
1349       ierr = SetTrajTLR(ts,s,stepnum,time,X);CHKERRQ(ierr);
1350       break;
1351     case REVOLVE_OFFLINE:
1352       ierr = SetTrajROF(ts,s,stepnum,time,X);CHKERRQ(ierr);
1353       break;
1354     case REVOLVE_ONLINE:
1355       ierr = SetTrajRON(ts,s,stepnum,time,X);CHKERRQ(ierr);
1356       break;
1357     case REVOLVE_MULTISTAGE:
1358       ierr = SetTrajRMS(ts,s,stepnum,time,X);CHKERRQ(ierr);
1359       break;
1360     default:
1361       break;
1362   }
1363   PetscFunctionReturn(0);
1364 }
1365 
1366 #undef __FUNCT__
1367 #define __FUNCT__ "TSTrajectoryGet_Memory"
1368 PetscErrorCode TSTrajectoryGet_Memory(TSTrajectory tj,TS ts,PetscInt stepnum,PetscReal *t)
1369 {
1370   Stack          *s = (Stack*)tj->data;
1371   PetscErrorCode ierr;
1372 
1373   PetscFunctionBegin;
1374   ierr = TSGetTotalSteps(ts,&stepnum);CHKERRQ(ierr);
1375   if (stepnum == 0) PetscFunctionReturn(0);
1376   switch (s->stype) {
1377     case NONE:
1378       ierr = GetTrajN(ts,s,stepnum);CHKERRQ(ierr);
1379       break;
1380     case TWO_LEVEL_NOREVOLVE:
1381       ierr = GetTrajTLNR(ts,s,stepnum);CHKERRQ(ierr);
1382       break;
1383     case TWO_LEVEL_REVOLVE:
1384       ierr = GetTrajTLR(ts,s,stepnum);CHKERRQ(ierr);
1385       break;
1386     case REVOLVE_OFFLINE:
1387       ierr = GetTrajROF(ts,s,stepnum);CHKERRQ(ierr);
1388       break;
1389     case REVOLVE_ONLINE:
1390       ierr = GetTrajRON(ts,s,stepnum);CHKERRQ(ierr);
1391       break;
1392     case REVOLVE_MULTISTAGE:
1393       ierr = GetTrajRMS(ts,s,stepnum);CHKERRQ(ierr);
1394       break;
1395     default:
1396       break;
1397   }
1398   PetscFunctionReturn(0);
1399 }
1400 
1401 #undef __FUNCT__
1402 #define __FUNCT__ "TSTrajectoryDestroy_Memory"
1403 PETSC_EXTERN PetscErrorCode TSTrajectoryDestroy_Memory(TSTrajectory tj)
1404 {
1405   Stack          *s = (Stack*)tj->data;
1406   PetscErrorCode ierr;
1407 
1408   PetscFunctionBegin;
1409   if (s->stype > TWO_LEVEL_NOREVOLVE) {
1410 #ifdef PETSC_HAVE_REVOLVE
1411     revolve_reset();
1412 #endif
1413   }
1414   if (s->stype == REVOLVE_ONLINE) {
1415 #ifdef PETSC_HAVE_REVOLVE
1416     s->top = s->stacksize-1;
1417 #endif
1418   }
1419   ierr = StackDestroy(&s);CHKERRQ(ierr);
1420   PetscFunctionReturn(0);
1421 }
1422 
1423 /*MC
1424       TSTRAJECTORYMEMORY - Stores each solution of the ODE/ADE in memory
1425 
1426   Level: intermediate
1427 
1428 .seealso:  TSTrajectoryCreate(), TS, TSTrajectorySetType()
1429 
1430 M*/
1431 #undef __FUNCT__
1432 #define __FUNCT__ "TSTrajectoryCreate_Memory"
1433 PETSC_EXTERN PetscErrorCode TSTrajectoryCreate_Memory(TSTrajectory tj,TS ts)
1434 {
1435   Stack          *s;
1436   PetscErrorCode ierr;
1437 
1438   PetscFunctionBegin;
1439   tj->ops->set            = TSTrajectorySet_Memory;
1440   tj->ops->get            = TSTrajectoryGet_Memory;
1441   tj->ops->setup          = TSTrajectorySetUp_Memory;
1442   tj->ops->destroy        = TSTrajectoryDestroy_Memory;
1443   tj->ops->setfromoptions = TSTrajectorySetFromOptions_Memory;
1444 
1445   ierr = PetscCalloc1(1,&s);CHKERRQ(ierr);
1446   s->stype        = NONE;
1447   s->max_cps_ram  = -1; /* -1 indicates that it is not set */
1448   s->max_cps_disk = -1; /* -1 indicates that it is not set */
1449   s->stride       = 0; /* if not zero, two-level checkpointing will be used */
1450   s->use_online   = PETSC_FALSE;
1451   s->save_stack   = PETSC_TRUE;
1452   s->solution_only= PETSC_TRUE;
1453   tj->data        = s;
1454   PetscFunctionReturn(0);
1455 }
1456