Actual source code: trajmemory.c

petsc-3.9.4 2018-09-11
Report Typos and Errors
  1:  #include <petsc/private/tsimpl.h>
  2:  #include <petscsys.h>
  3: #if defined(PETSC_HAVE_REVOLVE)
  4: #include <revolve_c.h>
  5: #endif

  7: PetscLogEvent TSTrajectory_DiskWrite, TSTrajectory_DiskRead;

  9: typedef enum {NONE,TWO_LEVEL_NOREVOLVE,TWO_LEVEL_REVOLVE,TWO_LEVEL_TWO_REVOLVE,REVOLVE_OFFLINE,REVOLVE_ONLINE,REVOLVE_MULTISTAGE} SchedulerType;

 11: typedef struct _StackElement {
 12:   PetscInt  stepnum;
 13:   Vec       X;
 14:   Vec       *Y;
 15:   PetscReal time;
 16:   PetscReal timeprev; /* for no solution_only mode */
 17:   PetscReal timenext; /* for solution_only mode */
 18: } *StackElement;

 20: #if defined(PETSC_HAVE_REVOLVE)
 21: typedef struct _RevolveCTX {
 22:   PetscBool reverseonestep;
 23:   PetscInt  where;
 24:   PetscInt  snaps_in;
 25:   PetscInt  stepsleft;
 26:   PetscInt  check;
 27:   PetscInt  oldcapo;
 28:   PetscInt  capo;
 29:   PetscInt  fine;
 30:   PetscInt  info;
 31: } RevolveCTX;
 32: #endif

 34: typedef struct _Stack {
 35:   PetscInt      stacksize;
 36:   PetscInt      top;
 37:   StackElement  *container;
 38:   PetscInt      numY;
 39:   PetscBool     solution_only;
 40:   PetscBool     use_dram;
 41: } Stack;

 43: typedef struct _DiskStack {
 44:   PetscInt  stacksize;
 45:   PetscInt  top;
 46:   PetscInt  *container;
 47: } DiskStack;

 49: typedef struct _TJScheduler {
 50:   SchedulerType stype;
 51: #if defined(PETSC_HAVE_REVOLVE)
 52:   RevolveCTX    *rctx,*rctx2;
 53:   PetscBool     use_online;
 54:   PetscInt      store_stride;
 55: #endif
 56:   PetscBool     recompute;
 57:   PetscBool     skip_trajectory;
 58:   PetscBool     save_stack;
 59:   PetscInt      max_cps_ram;  /* maximum checkpoints in RAM */
 60:   PetscInt      max_cps_disk; /* maximum checkpoints on disk */
 61:   PetscInt      stride;
 62:   PetscInt      total_steps;  /* total number of steps */
 63:   Stack         stack;
 64:   DiskStack     diskstack;
 65: } TJScheduler;

 67: static PetscErrorCode TurnForwardWithStepsize(TS ts,PetscReal nextstepsize)
 68: {
 69:     PetscReal      stepsize;

 73:     /* reverse the direction */
 74:     TSGetTimeStep(ts,&stepsize);
 75:     stepsize = nextstepsize;
 76:     TSSetTimeStep(ts,stepsize);
 77:     return(0);
 78: }

 80: static PetscErrorCode TurnForward(TS ts)
 81: {
 82:   PetscReal      stepsize;

 86:   /* reverse the direction */
 87:   TSGetTimeStep(ts,&stepsize);
 88:   TSSetTimeStep(ts,-stepsize);
 89:   return(0);
 90: }

 92: static PetscErrorCode TurnBackward(TS ts)
 93: {
 94:   PetscReal      stepsize;

 98:   /* reverse the direction */
 99:   stepsize = ts->ptime_prev-ts->ptime;
100:   TSSetTimeStep(ts,stepsize);
101:   return(0);
102: }

104: static PetscErrorCode ElementCreate(TS ts,Stack *stack,StackElement *e)
105: {
106:   Vec            X;
107:   Vec            *Y;

111:   if (stack->use_dram) {
112:     PetscMallocSetDRAM();
113:   }
114:   PetscCalloc1(1,e);
115:   TSGetSolution(ts,&X);
116:   VecDuplicate(X,&(*e)->X);
117:   if (stack->numY > 0 && !stack->solution_only) {
118:     TSGetStages(ts,&stack->numY,&Y);
119:     VecDuplicateVecs(Y[0],stack->numY,&(*e)->Y);
120:   }
121:   if (stack->use_dram) {
122:     PetscMallocResetDRAM();
123:   }
124:   return(0);
125: }

127: static PetscErrorCode ElementSet(TS ts,Stack *stack,StackElement *e,PetscInt stepnum,PetscReal time,Vec X)
128: {
129:   Vec            *Y;
130:   PetscInt       i;
131:   PetscReal      timeprev;

135:   VecCopy(X,(*e)->X);
136:   if (stack->numY > 0 && !stack->solution_only) {
137:     TSGetStages(ts,&stack->numY,&Y);
138:     for (i=0;i<stack->numY;i++) {
139:       VecCopy(Y[i],(*e)->Y[i]);
140:     }
141:   }
142:   (*e)->stepnum = stepnum;
143:   (*e)->time    = time;
144:   /* for consistency */
145:   if (stepnum == 0) {
146:     (*e)->timeprev = (*e)->time - ts->time_step;
147:   } else {
148:     TSGetPrevTime(ts,&timeprev);
149:     (*e)->timeprev = timeprev;
150:   }
151:   return(0);
152: }

154: static PetscErrorCode ElementDestroy(Stack *stack,StackElement e)
155: {

159:   if (stack->use_dram) {
160:     PetscMallocSetDRAM();
161:   }
162:   VecDestroy(&e->X);
163:   if (stack->numY > 0 && !stack->solution_only) {
164:     VecDestroyVecs(stack->numY,&e->Y);
165:   }
166:   PetscFree(e);
167:   if (stack->use_dram) {
168:     PetscMallocResetDRAM();
169:   }
170:   return(0);
171: }

173: static PetscErrorCode StackResize(Stack *stack,PetscInt newsize)
174: {
175:   StackElement   *newcontainer;
176:   PetscInt       i;

180:   PetscMalloc1(newsize*sizeof(StackElement),&newcontainer);
181:   for (i=0;i<stack->stacksize;i++) {
182:     newcontainer[i] = stack->container[i];
183:   }
184:   PetscFree(stack->container);
185:   stack->container = newcontainer;
186:   stack->stacksize = newsize;
187:   return(0);
188: }

190: static PetscErrorCode StackPush(Stack *stack,StackElement e)
191: {
193:   if (stack->top+1 >= stack->stacksize) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MEMC,"Maximum stack size (%D) exceeded",stack->stacksize);
194:   stack->container[++stack->top] = e;
195:   return(0);
196: }

198: static PetscErrorCode StackPop(Stack *stack,StackElement *e)
199: {
201:   *e = NULL;
202:   if (stack->top == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEMC,"Empty stack");
203:   *e = stack->container[stack->top--];
204:   return(0);
205: }

207: static PetscErrorCode StackTop(Stack *stack,StackElement *e)
208: {
210:   *e = stack->container[stack->top];
211:   return(0);
212: }

214: static PetscErrorCode StackCreate(Stack *stack,PetscInt size,PetscInt ny)
215: {

219:   stack->top  = -1;
220:   stack->numY = ny;

222:   PetscMalloc1(size*sizeof(StackElement),&stack->container);
223:   return(0);
224: }

226: static PetscErrorCode StackDestroy(Stack *stack)
227: {
228:   PetscInt       i,n;
229:   StackElement   e;

233:   if (stack->top > -1) {
234:     n = stack->top+1; /* number of elements in the stack */
235:     for (i=0;i<n;i++) {
236:       StackPop(stack,&e);
237:       ElementDestroy(stack,e);
238:     }
239:   }
240:   PetscFree(stack->container);
241:   return(0);
242: }

244: #if defined(PETSC_HAVE_REVOLVE)
245: static PetscErrorCode StackFind(Stack *stack,StackElement *e,PetscInt index)
246: {
248:   *e = stack->container[index];
249:   return(0);
250: }
251: #endif

253: static PetscErrorCode OutputBIN(MPI_Comm comm,const char *filename,PetscViewer *viewer)
254: {

258:   PetscViewerCreate(comm,viewer);
259:   PetscViewerSetType(*viewer,PETSCVIEWERBINARY);
260:   PetscViewerFileSetMode(*viewer,FILE_MODE_WRITE);
261:   PetscViewerFileSetName(*viewer,filename);
262:   return(0);
263: }

265: static PetscErrorCode WriteToDisk(PetscInt stepnum,PetscReal time,PetscReal timeprev,Vec X,Vec *Y,PetscInt numY,PetscBool solution_only,PetscViewer viewer)
266: {
267:   PetscInt       i;

271:   PetscViewerBinaryWrite(viewer,&stepnum,1,PETSC_INT,PETSC_FALSE);
272:   VecView(X,viewer);
273:   for (i=0;!solution_only && i<numY;i++) {
274:     VecView(Y[i],viewer);
275:   }
276:   PetscViewerBinaryWrite(viewer,&time,1,PETSC_REAL,PETSC_FALSE);
277:   PetscViewerBinaryWrite(viewer,&timeprev,1,PETSC_REAL,PETSC_FALSE);
278:   return(0);
279: }

281: static PetscErrorCode ReadFromDisk(PetscInt *stepnum,PetscReal *time,PetscReal *timeprev,Vec X,Vec *Y,PetscInt numY,PetscBool solution_only,PetscViewer viewer)
282: {
283:   PetscInt       i;

287:   PetscViewerBinaryRead(viewer,stepnum,1,NULL,PETSC_INT);
288:   VecLoad(X,viewer);
289:   for (i=0;!solution_only && i<numY;i++) {
290:     VecLoad(Y[i],viewer);
291:   }
292:   PetscViewerBinaryRead(viewer,time,1,NULL,PETSC_REAL);
293:   PetscViewerBinaryRead(viewer,timeprev,1,NULL,PETSC_REAL);
294:   return(0);
295: }

297: static PetscErrorCode StackDumpAll(TSTrajectory tj,TS ts,Stack *stack,PetscInt id)
298: {
299:   Vec            *Y;
300:   PetscInt       i;
301:   StackElement   e = NULL;
302:   PetscViewer    viewer;
303:   char           filename[PETSC_MAX_PATH_LEN];
305:   MPI_Comm       comm;

308:   PetscObjectGetComm((PetscObject)ts,&comm);
309:   if (tj->monitor) {
310:     PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);
311:     PetscViewerASCIIPrintf(tj->monitor,"Dump stack to file\n");
312:     PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);
313:   }
314:   if (id == 1) {
315:     PetscMPIInt rank;
316:     MPI_Comm_rank(comm,&rank);
317:     if (!rank) {
318:       PetscRMTree("SA-data");
319:       PetscMkdir("SA-data");
320:     }
321:   }
322:   PetscSNPrintf(filename,sizeof(filename),"SA-data/SA-STACK%06d.bin",id);
323:   OutputBIN(comm,filename,&viewer);
324:   for (i=0;i<stack->stacksize;i++) {
325:     e = stack->container[i];
326:     PetscLogEventBegin(TSTrajectory_DiskWrite,ts,0,0,0);
327:     WriteToDisk(e->stepnum,e->time,e->timeprev,e->X,e->Y,stack->numY,stack->solution_only,viewer);
328:     PetscLogEventEnd(TSTrajectory_DiskWrite,ts,0,0,0);
329:     ts->trajectory->diskwrites++;
330:   }
331:   /* 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 */
332:   TSGetStages(ts,&stack->numY,&Y);
333:   PetscLogEventBegin(TSTrajectory_DiskWrite,ts,0,0,0);
334:   WriteToDisk(ts->steps,ts->ptime,ts->ptime_prev,ts->vec_sol,Y,stack->numY,stack->solution_only,viewer);
335:   PetscLogEventEnd(TSTrajectory_DiskWrite,ts,0,0,0);
336:   ts->trajectory->diskwrites++;
337:   for (i=0;i<stack->stacksize;i++) {
338:     StackPop(stack,&e);
339:     ElementDestroy(stack,e);
340:   }
341:   PetscViewerDestroy(&viewer);
342:   return(0);
343: }

345: static PetscErrorCode StackLoadAll(TSTrajectory tj,TS ts,Stack *stack,PetscInt id)
346: {
347:   Vec            *Y;
348:   PetscInt       i;
349:   StackElement   e;
350:   PetscViewer    viewer;
351:   char           filename[PETSC_MAX_PATH_LEN];

355:   if (tj->monitor) {
356:     PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);
357:     PetscViewerASCIIPrintf(tj->monitor,"Load stack from file\n");
358:     PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);
359:   }
360:   PetscSNPrintf(filename,sizeof filename,"SA-data/SA-STACK%06d.bin",id);
361:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);
362:   for (i=0;i<stack->stacksize;i++) {
363:     ElementCreate(ts,stack,&e);
364:     StackPush(stack,e);
365:     PetscLogEventBegin(TSTrajectory_DiskRead,ts,0,0,0);
366:     ReadFromDisk(&e->stepnum,&e->time,&e->timeprev,e->X,e->Y,stack->numY,stack->solution_only,viewer);
367:     PetscLogEventEnd(TSTrajectory_DiskRead,ts,0,0,0);
368:     ts->trajectory->diskreads++;
369:   }
370:   /* load the last step into TS */
371:   TSGetStages(ts,&stack->numY,&Y);
372:   PetscLogEventBegin(TSTrajectory_DiskRead,ts,0,0,0);
373:   ReadFromDisk(&ts->steps,&ts->ptime,&ts->ptime_prev,ts->vec_sol,Y,stack->numY,stack->solution_only,viewer);
374:   PetscLogEventEnd(TSTrajectory_DiskRead,ts,0,0,0);
375:   ts->trajectory->diskreads++;
376:   TurnBackward(ts);
377:   PetscViewerDestroy(&viewer);
378:   return(0);
379: }

381: #if defined(PETSC_HAVE_REVOLVE)
382: static PetscErrorCode StackLoadLast(TSTrajectory tj,TS ts,Stack *stack,PetscInt id)
383: {
384:   Vec            *Y;
385:   PetscInt       size;
386:   PetscViewer    viewer;
387:   char           filename[PETSC_MAX_PATH_LEN];
388: #if defined(PETSC_HAVE_MPIIO)
389:   PetscBool      usempiio;
390: #endif
391:   int            fd;
392:   off_t          off,offset;

396:   if (tj->monitor) {
397:     PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);
398:     PetscViewerASCIIPrintf(tj->monitor,"Load last stack element from file\n");
399:     PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);
400:   }
401:   TSGetStages(ts,&stack->numY,&Y);
402:   VecGetSize(Y[0],&size);
403:   /* VecView writes to file two extra int's for class id and number of rows */
404:   off  = -((stack->solution_only?0:stack->numY)+1)*(size*PETSC_BINARY_SCALAR_SIZE+2*PETSC_BINARY_INT_SIZE)-PETSC_BINARY_INT_SIZE-2*PETSC_BINARY_SCALAR_SIZE;

406:   PetscSNPrintf(filename,sizeof filename,"SA-data/SA-STACK%06d.bin",id);
407:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);
408: #if defined(PETSC_HAVE_MPIIO)
409:   PetscViewerBinaryGetUseMPIIO(viewer,&usempiio);
410:   if (usempiio) {
411:     PetscViewerBinaryGetMPIIODescriptor(viewer,(MPI_File*)&fd);
412:     PetscBinarySynchronizedSeek(PETSC_COMM_WORLD,fd,off,PETSC_BINARY_SEEK_END,&offset);
413:   } else {
414: #endif
415:     PetscViewerBinaryGetDescriptor(viewer,&fd);
416:     PetscBinarySeek(fd,off,PETSC_BINARY_SEEK_END,&offset);
417: #if defined(PETSC_HAVE_MPIIO)
418:   }
419: #endif
420:   /* load the last step into TS */
421:   PetscLogEventBegin(TSTrajectory_DiskRead,ts,0,0,0);
422:   ReadFromDisk(&ts->steps,&ts->ptime,&ts->ptime_prev,ts->vec_sol,Y,stack->numY,stack->solution_only,viewer);
423:   PetscLogEventEnd(TSTrajectory_DiskRead,ts,0,0,0);
424:   ts->trajectory->diskreads++;
425:   TurnBackward(ts);
426:   PetscViewerDestroy(&viewer);
427:   return(0);
428: }
429: #endif

431: static PetscErrorCode DumpSingle(TSTrajectory tj,TS ts,Stack *stack,PetscInt id)
432: {
433:   Vec            *Y;
434:   PetscInt       stepnum;
435:   PetscViewer    viewer;
436:   char           filename[PETSC_MAX_PATH_LEN];
438:   MPI_Comm       comm;

441:   PetscObjectGetComm((PetscObject)ts,&comm);
442:   if (tj->monitor) {
443:     PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);
444:     PetscViewerASCIIPrintf(tj->monitor,"Load a single point from file\n");
445:     PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);
446:   }
447:   TSGetStepNumber(ts,&stepnum);
448:   if (id == 1) {
449:     PetscMPIInt rank;
450:     MPI_Comm_rank(comm,&rank);
451:     if (!rank) {
452:       PetscRMTree("SA-data");
453:       PetscMkdir("SA-data");
454:     }
455:   }
456:   PetscSNPrintf(filename,sizeof(filename),"SA-data/SA-CPS%06d.bin",id);
457:   OutputBIN(comm,filename,&viewer);

459:   TSGetStages(ts,&stack->numY,&Y);
460:   PetscLogEventBegin(TSTrajectory_DiskWrite,ts,0,0,0);
461:   WriteToDisk(stepnum,ts->ptime,ts->ptime_prev,ts->vec_sol,Y,stack->numY,stack->solution_only,viewer);
462:   PetscLogEventEnd(TSTrajectory_DiskWrite,ts,0,0,0);
463:   ts->trajectory->diskwrites++;

465:   PetscViewerDestroy(&viewer);
466:   return(0);
467: }

469: static PetscErrorCode LoadSingle(TSTrajectory tj,TS ts,Stack *stack,PetscInt id)
470: {
471:   Vec            *Y;
472:   PetscViewer    viewer;
473:   char           filename[PETSC_MAX_PATH_LEN];

477:   if (tj->monitor) {
478:     PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);
479:     PetscViewerASCIIPrintf(tj->monitor,"Load a single point from file\n");
480:     PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);
481:   }
482:   PetscSNPrintf(filename,sizeof filename,"SA-data/SA-CPS%06d.bin",id);
483:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);

485:   TSGetStages(ts,&stack->numY,&Y);
486:   PetscLogEventBegin(TSTrajectory_DiskRead,ts,0,0,0);
487:   ReadFromDisk(&ts->steps,&ts->ptime,&ts->ptime_prev,ts->vec_sol,Y,stack->numY,stack->solution_only,viewer);
488:   PetscLogEventEnd(TSTrajectory_DiskRead,ts,0,0,0);
489:   ts->trajectory->diskreads++;

491:   PetscViewerDestroy(&viewer);
492:   return(0);
493: }

495: static PetscErrorCode UpdateTS(TS ts,Stack *stack,StackElement e)
496: {
497:   Vec            *Y;
498:   PetscInt       i;

502:   VecCopy(e->X,ts->vec_sol);
503:   if (!stack->solution_only) {
504:     TSGetStages(ts,&stack->numY,&Y);
505:     for (i=0;i<stack->numY;i++) {
506:       VecCopy(e->Y[i],Y[i]);
507:     }
508:   }
509:   TSSetTimeStep(ts,e->timeprev-e->time); /* stepsize will be negative */
510:   ts->ptime      = e->time;
511:   ts->ptime_prev = e->timeprev;
512:   return(0);
513: }

515: static PetscErrorCode ReCompute(TS ts,TJScheduler *tjsch,PetscInt stepnumbegin,PetscInt stepnumend)
516: {
517:   Stack          *stack = &tjsch->stack;
518:   PetscInt       i;

522:   ts->steps = stepnumbegin; /* global step number */
523:   for (i=stepnumbegin;i<stepnumend;i++) { /* assume fixed step size */
524:     if (stack->solution_only && !tjsch->skip_trajectory) { /* revolve online need this */
525:       TSTrajectorySet(ts->trajectory,ts,ts->steps,ts->ptime,ts->vec_sol);
526:     }
527:     TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);
528:     TSStep(ts);
529:     if (!stack->solution_only && !tjsch->skip_trajectory) {
530:       TSTrajectorySet(ts->trajectory,ts,ts->steps,ts->ptime,ts->vec_sol);
531:     }
532:     TSEventHandler(ts);
533:     if (!ts->steprollback) {
534:       TSPostStep(ts);
535:     }
536:   }
537:   TurnBackward(ts);
538:   ts->trajectory->recomps += stepnumend-stepnumbegin; /* recomputation counter */
539:   ts->steps = stepnumend;
540:   return(0);
541: }

543: static PetscErrorCode TopLevelStore(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscInt localstepnum,PetscInt laststridesize,PetscBool *done)
544: {
545:   Stack          *stack = &tjsch->stack;
546:   DiskStack      *diskstack = &tjsch->diskstack;
547:   PetscInt       stridenum;

551:   *done = PETSC_FALSE;
552:   stridenum    = stepnum/tjsch->stride;
553:   /* make sure saved checkpoint id starts from 1
554:      skip last stride when using stridenum+1
555:      skip first stride when using stridenum */
556:   if (stack->solution_only) {
557:     if (tjsch->save_stack) {
558:       if (localstepnum == tjsch->stride-1 && stepnum < tjsch->total_steps-laststridesize) { /* current step will be saved without going through stack */
559:         StackDumpAll(tj,ts,stack,stridenum+1);
560:         if (tjsch->stype == TWO_LEVEL_TWO_REVOLVE) diskstack->container[++diskstack->top] = stridenum+1;
561:         *done = PETSC_TRUE;
562:       }
563:     } else {
564:       if (localstepnum == 0 && stepnum < tjsch->total_steps-laststridesize) {
565:         DumpSingle(tj,ts,stack,stridenum+1);
566:         if (tjsch->stype == TWO_LEVEL_TWO_REVOLVE) diskstack->container[++diskstack->top] = stridenum+1;
567:         *done = PETSC_TRUE;
568:       }
569:     }
570:   } else {
571:     if (tjsch->save_stack) {
572:       if (localstepnum == 0 && stepnum < tjsch->total_steps && stepnum != 0) { /* skip the first stride */
573:         StackDumpAll(tj,ts,stack,stridenum);
574:         if (tjsch->stype == TWO_LEVEL_TWO_REVOLVE) diskstack->container[++diskstack->top] = stridenum;
575:         *done = PETSC_TRUE;
576:       }
577:     } else {
578:       if (localstepnum == 1 && stepnum < tjsch->total_steps-laststridesize) {
579:         DumpSingle(tj,ts,stack,stridenum+1);
580:         if (tjsch->stype == TWO_LEVEL_TWO_REVOLVE) diskstack->container[++diskstack->top] = stridenum+1;
581:         *done = PETSC_TRUE;
582:       }
583:     }
584:   }
585:   return(0);
586: }

588: static PetscErrorCode SetTrajN(TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X)
589: {
590:   Stack          *stack = &tjsch->stack;
591:   StackElement   e;

595:   /* skip the last step */
596:   if (ts->reason) { /* only affect the forward run */
597:     /* update total_steps in the end of forward run */
598:     if (stepnum != tjsch->total_steps) tjsch->total_steps = stepnum;
599:     if (stack->solution_only) {
600:       /* get rid of the solution at second last step */
601:       StackPop(stack,&e);
602:       ElementDestroy(stack,e);
603:     }
604:     return(0);
605:   }
606:   /*  do not save trajectory at the recompute stage for solution_only mode */
607:   if (stack->solution_only && tjsch->recompute) return(0);
608:   /* skip the first step for no_solution_only mode */
609:   if (!stack->solution_only && stepnum == 0) return(0);

611:   /* resize the stack */
612:   if (stack->top+1 == stack->stacksize) {
613:     StackResize(stack,2*stack->stacksize);
614:   }
615:   /* update timenext for the previous step; necessary for step adaptivity */
616:   if (stack->top > -1) {
617:     StackTop(stack,&e);
618:     e->timenext = ts->ptime;
619:   }
620:   if (stepnum < stack->top) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_MEMC,"Illegal modification of a non-top stack element");
621:   ElementCreate(ts,stack,&e);
622:   ElementSet(ts,stack,&e,stepnum,time,X);
623:   StackPush(stack,e);
624:   return(0);
625: }

627: static PetscErrorCode GetTrajN(TS ts,TJScheduler *tjsch,PetscInt stepnum)
628: {
629:   Stack          *stack = &tjsch->stack;
630:   StackElement   e;

634:   if (stepnum == tjsch->total_steps) {
635:     TurnBackward(ts);
636:     return(0);
637:   }
638:   /* restore a checkpoint */
639:   StackTop(stack,&e);
640:   UpdateTS(ts,stack,e);
641:   if (stack->solution_only) {/* recompute one step */
642:     tjsch->recompute = PETSC_TRUE;
643:     TurnForwardWithStepsize(ts,e->timenext-e->time);
644:     ReCompute(ts,tjsch,e->stepnum,stepnum);
645:   }
646:   StackPop(stack,&e);
647:   ElementDestroy(stack,e);
648:   return(0);
649: }

651: static PetscErrorCode SetTrajTLNR(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X)
652: {
653:   Stack          *stack = &tjsch->stack;
654:   PetscInt       localstepnum,laststridesize;
655:   StackElement   e;
656:   PetscBool      done;

660:   if (!stack->solution_only && stepnum == 0) return(0);
661:   if (stack->solution_only && stepnum == tjsch->total_steps) return(0);
662:   if (tjsch->save_stack && tjsch->recompute) return(0);

664:   localstepnum = stepnum%tjsch->stride;
665:   /* (stride size-1) checkpoints are saved in each stride; an extra point is added by StackDumpAll() */
666:   laststridesize = tjsch->total_steps%tjsch->stride;
667:   if (!laststridesize) laststridesize = tjsch->stride;

669:   if (!tjsch->recompute) {
670:     TopLevelStore(tj,ts,tjsch,stepnum,localstepnum,laststridesize,&done);
671:     if (!tjsch->save_stack && stepnum < tjsch->total_steps-laststridesize) return(0);
672:   }
673:   if (!stack->solution_only && localstepnum == 0) return(0); /* skip last point in each stride at recompute stage or last stride */
674:   if (stack->solution_only && localstepnum == tjsch->stride-1) return(0); /* skip last step in each stride at recompute stage or last stride */

676:   ElementCreate(ts,stack,&e);
677:   ElementSet(ts,stack,&e,stepnum,time,X);
678:   StackPush(stack,e);
679:   return(0);
680: }

682: static PetscErrorCode GetTrajTLNR(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum)
683: {
684:   Stack          *stack = &tjsch->stack;
685:   PetscInt       id,localstepnum,laststridesize;
686:   StackElement   e;

690:   if (stepnum == tjsch->total_steps) {
691:     TurnBackward(ts);
692:     return(0);
693:   }

695:   localstepnum = stepnum%tjsch->stride;
696:   laststridesize = tjsch->total_steps%tjsch->stride;
697:   if (!laststridesize) laststridesize = tjsch->stride;
698:   if (stack->solution_only) {
699:     /* fill stack with info */
700:     if (localstepnum == 0 && tjsch->total_steps-stepnum >= laststridesize) {
701:       id = stepnum/tjsch->stride;
702:       if (tjsch->save_stack) {
703:         StackLoadAll(tj,ts,stack,id);
704:         tjsch->recompute = PETSC_TRUE;
705:         tjsch->skip_trajectory = PETSC_TRUE;
706:         TurnForward(ts);
707:         ReCompute(ts,tjsch,id*tjsch->stride-1,id*tjsch->stride);
708:         tjsch->skip_trajectory = PETSC_FALSE;
709:       } else {
710:         LoadSingle(tj,ts,stack,id);
711:         tjsch->recompute = PETSC_TRUE;
712:         TurnForward(ts);
713:         ReCompute(ts,tjsch,(id-1)*tjsch->stride,id*tjsch->stride);
714:       }
715:       return(0);
716:     }
717:     /* restore a checkpoint */
718:     StackPop(stack,&e);
719:     UpdateTS(ts,stack,e);
720:     tjsch->recompute = PETSC_TRUE;
721:     tjsch->skip_trajectory = PETSC_TRUE;
722:     TurnForward(ts);
723:     ReCompute(ts,tjsch,e->stepnum,stepnum);
724:     tjsch->skip_trajectory = PETSC_FALSE;
725:     ElementDestroy(stack,e);
726:   } else {
727:     /* fill stack with info */
728:     if (localstepnum == 0 && tjsch->total_steps-stepnum >= laststridesize) {
729:       id = stepnum/tjsch->stride;
730:       if (tjsch->save_stack) {
731:         StackLoadAll(tj,ts,stack,id);
732:       } else {
733:         LoadSingle(tj,ts,stack,id);
734:         ElementCreate(ts,stack,&e);
735:         ElementSet(ts,stack,&e,(id-1)*tjsch->stride+1,ts->ptime,ts->vec_sol);
736:         StackPush(stack,e);
737:         tjsch->recompute = PETSC_TRUE;
738:         TurnForward(ts);
739:         ReCompute(ts,tjsch,e->stepnum,id*tjsch->stride);
740:       }
741:       return(0);
742:     }
743:     /* restore a checkpoint */
744:     StackPop(stack,&e);
745:     UpdateTS(ts,stack,e);
746:     ElementDestroy(stack,e);
747:   }
748:   return(0);
749: }

751: #if defined(PETSC_HAVE_REVOLVE)
752: static PetscErrorCode printwhattodo(PetscViewer viewer,PetscInt whattodo,RevolveCTX *rctx,PetscInt shift)
753: {

757:   if (!viewer) return(0);

759:   switch(whattodo) {
760:     case 1:
761:       PetscViewerASCIIPrintf(viewer,"Advance from %D to %D\n",rctx->oldcapo+shift,rctx->capo+shift);
762:       break;
763:     case 2:
764:       PetscViewerASCIIPrintf(viewer,"Store in checkpoint number %D (located in RAM)\n",rctx->check);
765:       break;
766:     case 3:
767:       PetscViewerASCIIPrintf(viewer,"First turn: Initialize adjoints and reverse first step\n");
768:       break;
769:     case 4:
770:       PetscViewerASCIIPrintf(viewer,"Forward and reverse one step\n");
771:       break;
772:     case 5:
773:       PetscViewerASCIIPrintf(viewer,"Restore in checkpoint number %D (located in RAM)\n",rctx->check);
774:       break;
775:     case 7:
776:       PetscViewerASCIIPrintf(viewer,"Store in checkpoint number %D (located on disk)\n",rctx->check);
777:       break;
778:     case 8:
779:       PetscViewerASCIIPrintf(viewer,"Restore in checkpoint number %D (located on disk)\n",rctx->check);
780:       break;
781:     case -1:
782:       PetscViewerASCIIPrintf(viewer,"Error!");
783:       break;
784:   }
785:   return(0);
786: }

788: static PetscErrorCode printwhattodo2(PetscViewer viewer,PetscInt whattodo,RevolveCTX *rctx,PetscInt shift)
789: {

793:   if (!viewer) return(0);

795:   switch(whattodo) {
796:     case 1:
797:       PetscViewerASCIIPrintf(viewer,"[Top Level] Advance from stride %D to stride %D\n",rctx->oldcapo+shift,rctx->capo+shift);
798:       break;
799:     case 2:
800:       PetscViewerASCIIPrintf(viewer,"[Top Level] Store in checkpoint number %D\n",rctx->check);
801:       break;
802:     case 3:
803:       PetscViewerASCIIPrintf(viewer,"[Top Level] First turn: Initialize adjoints and reverse first stride\n");
804:       break;
805:     case 4:
806:       PetscViewerASCIIPrintf(viewer,"[Top Level] Forward and reverse one stride\n");
807:       break;
808:     case 5:
809:       PetscViewerASCIIPrintf(viewer,"[Top Level] Restore in checkpoint number %D\n",rctx->check);
810:       break;
811:     case 7:
812:       PetscViewerASCIIPrintf(viewer,"[Top Level] Store in top-level checkpoint number %D\n",rctx->check);
813:       break;
814:     case 8:
815:       PetscViewerASCIIPrintf(viewer,"[Top Level] Restore in top-level checkpoint number %D\n",rctx->check);
816:       break;
817:     case -1:
818:       PetscViewerASCIIPrintf(viewer,"[Top Level] Error!");
819:       break;
820:   }
821:   return(0);
822: }

824: static PetscErrorCode InitRevolve(PetscInt fine,PetscInt snaps,RevolveCTX *rctx)
825: {
827:   revolve_reset();
828:   revolve_create_offline(fine,snaps);
829:   rctx->snaps_in       = snaps;
830:   rctx->fine           = fine;
831:   rctx->check          = 0;
832:   rctx->capo           = 0;
833:   rctx->reverseonestep = PETSC_FALSE;
834:   /* check stepsleft? */
835:   return(0);
836: }

838: static PetscErrorCode FastForwardRevolve(RevolveCTX *rctx)
839: {
840:   PetscInt whattodo;

843:   whattodo = 0;
844:   while(whattodo!=3) { /* we have to fast forward revolve to the beginning of the backward sweep due to unfriendly revolve interface */
845:     whattodo = revolve_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where);
846:   }
847:   return(0);
848: }

850: static PetscErrorCode ApplyRevolve(PetscViewer viewer,SchedulerType stype,RevolveCTX *rctx,PetscInt total_steps,PetscInt stepnum,PetscInt localstepnum,PetscBool toplevel,PetscInt *store)
851: {
853:   PetscInt       shift,whattodo;

856:   *store = 0;
857:   if (rctx->stepsleft > 0) { /* advance the solution without checkpointing anything as Revolve requires */
858:     rctx->stepsleft--;
859:     return(0);
860:   }
861:   /* let Revolve determine what to do next */
862:   shift         = stepnum-localstepnum;
863:   rctx->oldcapo = rctx->capo;
864:   rctx->capo    = localstepnum;

866:   if (!toplevel) whattodo = revolve_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where);
867:   else whattodo = revolve2_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where);
868:   if (stype == REVOLVE_ONLINE && whattodo == 8) whattodo = 5;
869:   if (stype == REVOLVE_ONLINE && whattodo == 7) whattodo = 2;
870:   if (!toplevel) {printwhattodo(viewer,whattodo,rctx,shift);}
871:   else {printwhattodo2(viewer,whattodo,rctx,shift);}
872:   if (whattodo == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in the Revolve library");
873:   if (whattodo == 1) { /* advance some time steps */
874:     if (stype == REVOLVE_ONLINE && rctx->capo >= total_steps-1) {
875:       revolve_turn(total_steps,&rctx->capo,&rctx->fine);
876:       if (!toplevel) {printwhattodo(viewer,whattodo,rctx,shift);}
877:       else {printwhattodo2(viewer,whattodo,rctx,shift);}
878:     }
879:     rctx->stepsleft = rctx->capo-rctx->oldcapo-1;
880:   }
881:   if (whattodo == 3 || whattodo == 4) { /* ready for a reverse step */
882:     rctx->reverseonestep = PETSC_TRUE;
883:   }
884:   if (whattodo == 5) { /* restore a checkpoint and ask Revolve what to do next */
885:     rctx->oldcapo = rctx->capo;
886:     if (!toplevel) whattodo = revolve_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where); /* must return 1 or 3 or 4*/
887:     else whattodo = revolve2_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where);
888:     if (!toplevel) {printwhattodo(viewer,whattodo,rctx,shift);}
889:     else {printwhattodo2(viewer,whattodo,rctx,shift);}
890:     if (whattodo == 3 || whattodo == 4) rctx->reverseonestep = PETSC_TRUE;
891:     if (whattodo == 1) rctx->stepsleft = rctx->capo-rctx->oldcapo;
892:   }
893:   if (whattodo == 7) { /* save the checkpoint to disk */
894:     *store = 2;
895:     rctx->oldcapo = rctx->capo;
896:     whattodo = revolve_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where); /* must return 1 */
897:     printwhattodo(viewer,whattodo,rctx,shift);
898:     rctx->stepsleft = rctx->capo-rctx->oldcapo-1;
899:   }
900:   if (whattodo == 2) { /* store a checkpoint to RAM and ask Revolve how many time steps to advance next */
901:     *store = 1;
902:     rctx->oldcapo = rctx->capo;
903:     if (!toplevel) whattodo = revolve_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where); /* must return 1 */
904:     else whattodo = revolve2_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where);
905:     if (!toplevel) {printwhattodo(viewer,whattodo,rctx,shift);}
906:     else {printwhattodo2(viewer,whattodo,rctx,shift);}
907:     if (stype == REVOLVE_ONLINE && rctx->capo >= total_steps-1) {
908:       revolve_turn(total_steps,&rctx->capo,&rctx->fine);
909:       printwhattodo(viewer,whattodo,rctx,shift);
910:     }
911:     rctx->stepsleft = rctx->capo-rctx->oldcapo-1;
912:   }
913:   return(0);
914: }

916: static PetscErrorCode SetTrajROF(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X)
917: {
918:   Stack          *stack = &tjsch->stack;
919:   PetscInt       store;
920:   StackElement   e;

924:   if (!stack->solution_only && stepnum == 0) return(0);
925:   if (stack->solution_only && stepnum == tjsch->total_steps) return(0);
926:   ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,tjsch->total_steps,stepnum,stepnum,PETSC_FALSE,&store);
927:   if (store == 1) {
928:     if (stepnum < stack->top) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_MEMC,"Illegal modification of a non-top stack element");
929:     ElementCreate(ts,stack,&e);
930:     ElementSet(ts,stack,&e,stepnum,time,X);
931:     StackPush(stack,e);
932:   }
933:   return(0);
934: }

936: static PetscErrorCode GetTrajROF(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum)
937: {
938:   Stack          *stack = &tjsch->stack;
939:   PetscInt       whattodo,shift,store;
940:   StackElement   e;

944:   if (stepnum == 0 || stepnum == tjsch->total_steps) {
945:     TurnBackward(ts);
946:     tjsch->rctx->reverseonestep = PETSC_FALSE;
947:     return(0);
948:   }
949:   /* restore a checkpoint */
950:   StackTop(stack,&e);
951:   UpdateTS(ts,stack,e);
952:   if (stack->solution_only) { /* start with restoring a checkpoint */
953:     tjsch->rctx->capo = stepnum;
954:     tjsch->rctx->oldcapo = tjsch->rctx->capo;
955:     shift = 0;
956:     whattodo = revolve_action(&tjsch->rctx->check,&tjsch->rctx->capo,&tjsch->rctx->fine,tjsch->rctx->snaps_in,&tjsch->rctx->info,&tjsch->rctx->where);
957:     printwhattodo(tj->monitor,whattodo,tjsch->rctx,shift);
958:   } else { /* 2 revolve actions: restore a checkpoint and then advance */
959:     ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,tjsch->total_steps,stepnum,stepnum,PETSC_FALSE,&store);
960:     if (tj->monitor) {
961:       PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);
962:       PetscViewerASCIIPrintf(tj->monitor,"Skip the step from %D to %D (stage values already checkpointed)\n",tjsch->rctx->oldcapo,tjsch->rctx->oldcapo+1);
963:       PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);
964:     }
965:     if (!tjsch->rctx->reverseonestep && tjsch->rctx->stepsleft > 0) tjsch->rctx->stepsleft--;
966:   }
967:   if (stack->solution_only || (!stack->solution_only && e->stepnum < stepnum)) {
968:     tjsch->recompute = PETSC_TRUE;
969:     TurnForward(ts);
970:     ReCompute(ts,tjsch,e->stepnum,stepnum);
971:   }
972:   if ((stack->solution_only && e->stepnum+1 == stepnum) || (!stack->solution_only && e->stepnum == stepnum)) {
973:     StackPop(stack,&e);
974:     ElementDestroy(stack,e);
975:   }
976:   tjsch->rctx->reverseonestep = PETSC_FALSE;
977:   return(0);
978: }

980: static PetscErrorCode SetTrajRON(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X)
981: {
982:   Stack          *stack = &tjsch->stack;
983:   Vec            *Y;
984:   PetscInt       i,store;
985:   PetscReal      timeprev;
986:   StackElement   e;
987:   RevolveCTX     *rctx = tjsch->rctx;

991:   if (!stack->solution_only && stepnum == 0) return(0);
992:   if (stack->solution_only && stepnum == tjsch->total_steps) return(0);
993:   ApplyRevolve(tj->monitor,tjsch->stype,rctx,tjsch->total_steps,stepnum,stepnum,PETSC_FALSE,&store);
994:   if (store == 1) {
995:     if (rctx->check != stack->top+1) { /* overwrite some non-top checkpoint in the stack */
996:       StackFind(stack,&e,rctx->check);
997:       VecCopy(X,e->X);
998:       if (stack->numY > 0 && !stack->solution_only) {
999:         TSGetStages(ts,&stack->numY,&Y);
1000:         for (i=0;i<stack->numY;i++) {
1001:           VecCopy(Y[i],e->Y[i]);
1002:         }
1003:       }
1004:       e->stepnum  = stepnum;
1005:       e->time     = time;
1006:       TSGetPrevTime(ts,&timeprev);
1007:       e->timeprev = timeprev;
1008:     } else {
1009:       if (stepnum < stack->top) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_MEMC,"Illegal modification of a non-top stack element");
1010:       ElementCreate(ts,stack,&e);
1011:       ElementSet(ts,stack,&e,stepnum,time,X);
1012:       StackPush(stack,e);
1013:     }
1014:   }
1015:   return(0);
1016: }

1018: static PetscErrorCode GetTrajRON(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum)
1019: {
1020:   Stack          *stack = &tjsch->stack;
1021:   PetscInt       whattodo,shift;
1022:   StackElement   e;

1026:   if (stepnum == 0 || stepnum == tjsch->total_steps) {
1027:     TurnBackward(ts);
1028:     tjsch->rctx->reverseonestep = PETSC_FALSE;
1029:     return(0);
1030:   }
1031:   tjsch->rctx->capo = stepnum;
1032:   tjsch->rctx->oldcapo = tjsch->rctx->capo;
1033:   shift = 0;
1034:   whattodo = revolve_action(&tjsch->rctx->check,&tjsch->rctx->capo,&tjsch->rctx->fine,tjsch->rctx->snaps_in,&tjsch->rctx->info,&tjsch->rctx->where); /* whattodo=restore */
1035:   if (whattodo == 8) whattodo = 5;
1036:   printwhattodo(tj->monitor,whattodo,tjsch->rctx,shift);
1037:   /* restore a checkpoint */
1038:   StackFind(stack,&e,tjsch->rctx->check);
1039:   UpdateTS(ts,stack,e);
1040:   if (!stack->solution_only) { /* whattodo must be 5 */
1041:     /* ask Revolve what to do next */
1042:     tjsch->rctx->oldcapo = tjsch->rctx->capo;
1043:     whattodo = revolve_action(&tjsch->rctx->check,&tjsch->rctx->capo,&tjsch->rctx->fine,tjsch->rctx->snaps_in,&tjsch->rctx->info,&tjsch->rctx->where); /* must return 1 or 3 or 4*/
1044:     printwhattodo(tj->monitor,whattodo,tjsch->rctx,shift);
1045:     if (whattodo == 3 || whattodo == 4) tjsch->rctx->reverseonestep = PETSC_TRUE;
1046:     if (whattodo == 1) tjsch->rctx->stepsleft = tjsch->rctx->capo-tjsch->rctx->oldcapo;
1047:     if (tj->monitor) {
1048:       PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);
1049:       PetscViewerASCIIPrintf(tj->monitor,"Skip the step from %D to %D (stage values already checkpointed)\n",tjsch->rctx->oldcapo,tjsch->rctx->oldcapo+1);
1050:       PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);
1051:     }
1052:     if (!tjsch->rctx->reverseonestep && tjsch->rctx->stepsleft > 0) tjsch->rctx->stepsleft--;
1053:   }
1054:   if (stack->solution_only || (!stack->solution_only && e->stepnum < stepnum)) {
1055:     tjsch->recompute = PETSC_TRUE;
1056:     TurnForward(ts);
1057:     ReCompute(ts,tjsch,e->stepnum,stepnum);
1058:   }
1059:   tjsch->rctx->reverseonestep = PETSC_FALSE;
1060:   return(0);
1061: }

1063: static PetscErrorCode SetTrajTLR(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X)
1064: {
1065:   Stack          *stack = &tjsch->stack;
1066:   PetscInt       store,localstepnum,laststridesize;
1067:   StackElement   e;
1068:   PetscBool      done = PETSC_FALSE;

1072:   if (!stack->solution_only && stepnum == 0) return(0);
1073:   if (stack->solution_only && stepnum == tjsch->total_steps) return(0);

1075:   localstepnum = stepnum%tjsch->stride;
1076:   laststridesize = tjsch->total_steps%tjsch->stride;
1077:   if (!laststridesize) laststridesize = tjsch->stride;

1079:   if (!tjsch->recompute) {
1080:     TopLevelStore(tj,ts,tjsch,stepnum,localstepnum,laststridesize,&done);
1081:     /* revolve is needed for the last stride; different starting points for last stride between solutin_only and !solutin_only */
1082:     if (!stack->solution_only && !tjsch->save_stack && stepnum <= tjsch->total_steps-laststridesize) return(0);
1083:     if (stack->solution_only && !tjsch->save_stack && stepnum < tjsch->total_steps-laststridesize) return(0);
1084:   }
1085:   if (tjsch->save_stack && done) {
1086:     InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);
1087:     return(0);
1088:   }
1089:   if (laststridesize < tjsch->stride) {
1090:     if (stack->solution_only && stepnum == tjsch->total_steps-laststridesize && !tjsch->recompute) { /* step tjsch->total_steps-laststridesize-1 is skipped, but the next step is not */
1091:       InitRevolve(laststridesize,tjsch->max_cps_ram,tjsch->rctx);
1092:     }
1093:     if (!stack->solution_only && stepnum == tjsch->total_steps-laststridesize+1 && !tjsch->recompute) { /* step tjsch->total_steps-laststridesize is skipped, but the next step is not */
1094:       InitRevolve(laststridesize,tjsch->max_cps_ram,tjsch->rctx);
1095:     }
1096:   }
1097:   ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,tjsch->total_steps,stepnum,localstepnum,PETSC_FALSE,&store);
1098:   if (store == 1) {
1099:     if (localstepnum < stack->top) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_MEMC,"Illegal modification of a non-top stack element");
1100:     ElementCreate(ts,stack,&e);
1101:     ElementSet(ts,stack,&e,stepnum,time,X);
1102:     StackPush(stack,e);
1103:   }
1104:   return(0);
1105: }

1107: static PetscErrorCode GetTrajTLR(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum)
1108: {
1109:   Stack          *stack = &tjsch->stack;
1110:   PetscInt       whattodo,shift;
1111:   PetscInt       localstepnum,stridenum,laststridesize,store;
1112:   StackElement   e;

1116:   localstepnum = stepnum%tjsch->stride;
1117:   stridenum    = stepnum/tjsch->stride;
1118:   if (stepnum == tjsch->total_steps) {
1119:     TurnBackward(ts);
1120:     tjsch->rctx->reverseonestep = PETSC_FALSE;
1121:     return(0);
1122:   }
1123:   laststridesize = tjsch->total_steps%tjsch->stride;
1124:   if (!laststridesize) laststridesize = tjsch->stride;
1125:   if (stack->solution_only) {
1126:     /* fill stack */
1127:     if (localstepnum == 0 && stepnum <= tjsch->total_steps-laststridesize) {
1128:       if (tjsch->save_stack) {
1129:         StackLoadAll(tj,ts,stack,stridenum);
1130:         InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);
1131:         FastForwardRevolve(tjsch->rctx);
1132:         tjsch->recompute = PETSC_TRUE;
1133:         tjsch->skip_trajectory = PETSC_TRUE;
1134:         TurnForward(ts);
1135:         ReCompute(ts,tjsch,stridenum*tjsch->stride-1,stridenum*tjsch->stride);
1136:         tjsch->skip_trajectory = PETSC_FALSE;
1137:       } else {
1138:         LoadSingle(tj,ts,stack,stridenum);
1139:         InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);
1140:         tjsch->recompute = PETSC_TRUE;
1141:         TurnForward(ts);
1142:         ReCompute(ts,tjsch,(stridenum-1)*tjsch->stride,stridenum*tjsch->stride);
1143:       }
1144:       return(0);
1145:     }
1146:     /* restore a checkpoint */
1147:     StackTop(stack,&e);
1148:     UpdateTS(ts,stack,e);
1149:     /* start with restoring a checkpoint */
1150:     tjsch->rctx->capo = stepnum;
1151:     tjsch->rctx->oldcapo = tjsch->rctx->capo;
1152:     shift = stepnum-localstepnum;
1153:     whattodo = revolve_action(&tjsch->rctx->check,&tjsch->rctx->capo,&tjsch->rctx->fine,tjsch->rctx->snaps_in,&tjsch->rctx->info,&tjsch->rctx->where);
1154:     printwhattodo(tj->monitor,whattodo,tjsch->rctx,shift);
1155:     tjsch->recompute = PETSC_TRUE;
1156:     TurnForward(ts);
1157:     ReCompute(ts,tjsch,e->stepnum,stepnum);
1158:     if (e->stepnum+1 == stepnum) {
1159:       StackPop(stack,&e);
1160:       ElementDestroy(stack,e);
1161:     }
1162:   } else {
1163:     /* fill stack with info */
1164:     if (localstepnum == 0 && tjsch->total_steps-stepnum >= laststridesize) {
1165:       if (tjsch->save_stack) {
1166:         StackLoadAll(tj,ts,stack,stridenum);
1167:         InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);
1168:         FastForwardRevolve(tjsch->rctx);
1169:       } else {
1170:         LoadSingle(tj,ts,stack,stridenum);
1171:         InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);
1172:         ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,tjsch->total_steps,(stridenum-1)*tjsch->stride+1,1,PETSC_FALSE,&store);
1173:         if (tj->monitor) {
1174:           PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);
1175:           PetscViewerASCIIPrintf(tj->monitor,"Skip the step from %D to %D (stage values already checkpointed)\n",(stridenum-1)*tjsch->stride+tjsch->rctx->oldcapo,(stridenum-1)*tjsch->stride+tjsch->rctx->oldcapo+1);
1176:           PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);
1177:         }
1178:         ElementCreate(ts,stack,&e);
1179:         ElementSet(ts,stack,&e,(stridenum-1)*tjsch->stride+1,ts->ptime,ts->vec_sol);
1180:         StackPush(stack,e);
1181:         tjsch->recompute = PETSC_TRUE;
1182:         TurnForward(ts);
1183:         ReCompute(ts,tjsch,e->stepnum,stridenum*tjsch->stride);
1184:       }
1185:       return(0);
1186:     }
1187:     /* restore a checkpoint */
1188:     StackTop(stack,&e);
1189:     UpdateTS(ts,stack,e);
1190:     /* 2 revolve actions: restore a checkpoint and then advance */
1191:     ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,tjsch->total_steps,stepnum,localstepnum,PETSC_FALSE,&store);
1192:     if (tj->monitor) {
1193:       PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);
1194:       PetscViewerASCIIPrintf(tj->monitor,"Skip the step from %D to %D (stage values already checkpointed)\n",stepnum-localstepnum+tjsch->rctx->oldcapo,stepnum-localstepnum+tjsch->rctx->oldcapo+1);
1195:       PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);
1196:     }
1197:     if (!tjsch->rctx->reverseonestep && tjsch->rctx->stepsleft > 0) tjsch->rctx->stepsleft--;
1198:     if (e->stepnum < stepnum) {
1199:       tjsch->recompute = PETSC_TRUE;
1200:       TurnForward(ts);
1201:       ReCompute(ts,tjsch,e->stepnum,stepnum);
1202:     }
1203:     if (e->stepnum == stepnum) {
1204:       StackPop(stack,&e);
1205:       ElementDestroy(stack,e);
1206:     }
1207:   }
1208:   tjsch->rctx->reverseonestep = PETSC_FALSE;
1209:   return(0);
1210: }

1212: static PetscErrorCode SetTrajTLTR(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X)
1213: {
1214:   Stack          *stack = &tjsch->stack;
1215:   PetscInt       store,localstepnum,stridenum,laststridesize;
1216:   StackElement   e;
1217:   PetscBool      done = PETSC_FALSE;

1221:   if (!stack->solution_only && stepnum == 0) return(0);
1222:   if (stack->solution_only && stepnum == tjsch->total_steps) return(0);

1224:   localstepnum = stepnum%tjsch->stride; /* index at the bottom level (inside a stride) */
1225:   stridenum    = stepnum/tjsch->stride; /* index at the top level */
1226:   laststridesize = tjsch->total_steps%tjsch->stride;
1227:   if (!laststridesize) laststridesize = tjsch->stride;
1228:   if (stack->solution_only && localstepnum == 0 && !tjsch->rctx2->reverseonestep) {
1229:     ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx2,(tjsch->total_steps+tjsch->stride-1)/tjsch->stride,stridenum,stridenum,PETSC_TRUE,&tjsch->store_stride);
1230:     if (laststridesize < tjsch->stride && stepnum == tjsch->total_steps-laststridesize) {
1231:       InitRevolve(laststridesize,tjsch->max_cps_ram,tjsch->rctx);
1232:     }
1233:   }
1234:   if (!stack->solution_only && localstepnum == 1 && !tjsch->rctx2->reverseonestep) {
1235:     ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx2,(tjsch->total_steps+tjsch->stride-1)/tjsch->stride,stridenum,stridenum,PETSC_TRUE,&tjsch->store_stride);
1236:     if (laststridesize < tjsch->stride && stepnum == tjsch->total_steps-laststridesize+1) {
1237:       InitRevolve(laststridesize,tjsch->max_cps_ram,tjsch->rctx);
1238:     }
1239:   }
1240:   if (tjsch->store_stride) {
1241:     TopLevelStore(tj,ts,tjsch,stepnum,localstepnum,laststridesize,&done);
1242:     if (done) {
1243:       InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);
1244:       return(0);
1245:     }
1246:   }
1247:   if (stepnum < tjsch->total_steps-laststridesize) {
1248:     if (tjsch->save_stack && !tjsch->store_stride && !tjsch->rctx2->reverseonestep) return(0); /* store or forward-and-reverse at top level trigger revolve at bottom level */
1249:     if (!tjsch->save_stack && !tjsch->rctx2->reverseonestep) return(0); /* store operation does not require revolve be called at bottom level */
1250:   }
1251:   /* Skipping stepnum=0 for !stack->only is enough for TLR, but not for TLTR. Here we skip the first step for each stride so that the top-level revolve is applied (always at localstepnum=1) ahead of the bottom-level revolve */
1252:   if (!stack->solution_only && localstepnum == 0 && stepnum != tjsch->total_steps && !tjsch->recompute) return(0);
1253:   ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,tjsch->total_steps,stepnum,localstepnum,PETSC_FALSE,&store);
1254:   if (store == 1) {
1255:     if (localstepnum < stack->top) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_MEMC,"Illegal modification of a non-top stack element");
1256:     ElementCreate(ts,stack,&e);
1257:     ElementSet(ts,stack,&e,stepnum,time,X);
1258:     StackPush(stack,e);
1259:   }
1260:   return(0);
1261: }

1263: static PetscErrorCode GetTrajTLTR(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum)
1264: {
1265:   Stack          *stack = &tjsch->stack;
1266:   DiskStack      *diskstack = &tjsch->diskstack;
1267:   PetscInt       whattodo,shift;
1268:   PetscInt       localstepnum,stridenum,restoredstridenum,laststridesize,store;
1269:   StackElement   e;

1273:   localstepnum = stepnum%tjsch->stride;
1274:   stridenum    = stepnum/tjsch->stride;
1275:   if (stepnum == tjsch->total_steps) {
1276:     TurnBackward(ts);
1277:     tjsch->rctx->reverseonestep = PETSC_FALSE;
1278:     return(0);
1279:   }
1280:   laststridesize = tjsch->total_steps%tjsch->stride;
1281:   if (!laststridesize) laststridesize = tjsch->stride;
1282:   /*
1283:    Last stride can be adjoined directly. All the other strides require that the stack in memory be ready before an adjoint step is taken (at the end of each stride). The following two cases need to be addressed differently:
1284:      Case 1 (save_stack)
1285:        Restore a disk checkpoint; update TS with the last element in the restored data; recompute to the current point.
1286:      Case 2 (!save_stack)
1287:        Restore a disk checkpoint; update TS with the restored point; recompute to the current point.
1288:   */
1289:   if (localstepnum == 0 && stepnum <= tjsch->total_steps-laststridesize) {
1290:     /* restore the top element in the stack for disk checkpoints */
1291:     restoredstridenum = diskstack->container[diskstack->top];
1292:     tjsch->rctx2->reverseonestep = PETSC_FALSE;
1293:     /* top-level revolve must be applied before current step, just like the solution_only mode for single-level revolve */
1294:     if (!tjsch->save_stack && stack->solution_only) { /* start with restoring a checkpoint */
1295:       tjsch->rctx2->capo = stridenum;
1296:       tjsch->rctx2->oldcapo = tjsch->rctx2->capo;
1297:       shift = 0;
1298:       whattodo = revolve2_action(&tjsch->rctx2->check,&tjsch->rctx2->capo,&tjsch->rctx2->fine,tjsch->rctx2->snaps_in,&tjsch->rctx2->info,&tjsch->rctx2->where);
1299:       printwhattodo2(tj->monitor,whattodo,tjsch->rctx2,shift);
1300:     } else { /* 2 revolve actions: restore a checkpoint and then advance */
1301:       ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx2,(tjsch->total_steps+tjsch->stride-1)/tjsch->stride,stridenum,stridenum,PETSC_TRUE,&tjsch->store_stride);
1302:       if (tj->monitor) {
1303:         PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);
1304:         PetscViewerASCIIPrintf(tj->monitor,"[Top Level] Skip the stride from %D to %D (stage values already checkpointed)\n",tjsch->rctx2->oldcapo,tjsch->rctx2->oldcapo+1);
1305:         PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);
1306:       }
1307:       if (!tjsch->rctx2->reverseonestep && tjsch->rctx2->stepsleft > 0) tjsch->rctx2->stepsleft--;
1308:     }
1309:     /* fill stack */
1310:     if (stack->solution_only) {
1311:       if (tjsch->save_stack) {
1312:         if (restoredstridenum < stridenum) {
1313:           StackLoadLast(tj,ts,stack,restoredstridenum);
1314:         } else {
1315:           StackLoadAll(tj,ts,stack,restoredstridenum);
1316:         }
1317:         /* recompute one step ahead */
1318:         tjsch->recompute = PETSC_TRUE;
1319:         tjsch->skip_trajectory = PETSC_TRUE;
1320:         TurnForward(ts);
1321:         ReCompute(ts,tjsch,stridenum*tjsch->stride-1,stridenum*tjsch->stride);
1322:         tjsch->skip_trajectory = PETSC_FALSE;
1323:         if (restoredstridenum < stridenum) {
1324:           InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);
1325:           tjsch->recompute = PETSC_TRUE;
1326:           TurnForward(ts);
1327:           ReCompute(ts,tjsch,restoredstridenum*tjsch->stride,stepnum);
1328:         } else { /* stack ready, fast forward revolve status */
1329:           InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);
1330:           FastForwardRevolve(tjsch->rctx);
1331:         }
1332:       } else {
1333:         LoadSingle(tj,ts,stack,restoredstridenum);
1334:         InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);
1335:         tjsch->recompute = PETSC_TRUE;
1336:         TurnForward(ts);
1337:         ReCompute(ts,tjsch,(restoredstridenum-1)*tjsch->stride,stepnum);
1338:       }
1339:     } else {
1340:       if (tjsch->save_stack) {
1341:         if (restoredstridenum < stridenum) {
1342:           StackLoadLast(tj,ts,stack,restoredstridenum);
1343:           /* reset revolve */
1344:           InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);
1345:           tjsch->recompute = PETSC_TRUE;
1346:           TurnForward(ts);
1347:           ReCompute(ts,tjsch,restoredstridenum*tjsch->stride,stepnum);
1348:         } else { /* stack ready, fast forward revolve status */
1349:           StackLoadAll(tj,ts,stack,restoredstridenum);
1350:           InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);
1351:           FastForwardRevolve(tjsch->rctx);
1352:         }
1353:       } else {
1354:         LoadSingle(tj,ts,stack,restoredstridenum);
1355:         InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);
1356:         /* push first element to stack */
1357:         if (tjsch->store_stride || tjsch->rctx2->reverseonestep) {
1358:           shift = (restoredstridenum-1)*tjsch->stride-localstepnum;
1359:           ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,tjsch->total_steps,(restoredstridenum-1)*tjsch->stride+1,1,PETSC_FALSE,&store);
1360:           if (tj->monitor) {
1361:             PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);
1362:             PetscViewerASCIIPrintf(tj->monitor,"Skip the step from %D to %D (stage values already checkpointed)\n",(restoredstridenum-1)*tjsch->stride,(restoredstridenum-1)*tjsch->stride+1);
1363:             PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);
1364:           }
1365:           ElementCreate(ts,stack,&e);
1366:           ElementSet(ts,stack,&e,(restoredstridenum-1)*tjsch->stride+1,ts->ptime,ts->vec_sol);
1367:           StackPush(stack,e);
1368:         }
1369:         tjsch->recompute = PETSC_TRUE;
1370:         TurnForward(ts);
1371:         ReCompute(ts,tjsch,(restoredstridenum-1)*tjsch->stride+1,stepnum);
1372:       }
1373:     }
1374:     if (restoredstridenum == stridenum) diskstack->top--;
1375:     tjsch->rctx->reverseonestep = PETSC_FALSE;
1376:     return(0);
1377:   }

1379:   if (stack->solution_only) {
1380:     /* restore a checkpoint */
1381:     StackTop(stack,&e);
1382:     UpdateTS(ts,stack,e);
1383:     /* start with restoring a checkpoint */
1384:     tjsch->rctx->capo = stepnum;
1385:     tjsch->rctx->oldcapo = tjsch->rctx->capo;
1386:     shift = stepnum-localstepnum;
1387:     whattodo = revolve_action(&tjsch->rctx->check,&tjsch->rctx->capo,&tjsch->rctx->fine,tjsch->rctx->snaps_in,&tjsch->rctx->info,&tjsch->rctx->where);
1388:     printwhattodo(tj->monitor,whattodo,tjsch->rctx,shift);
1389:     tjsch->recompute = PETSC_TRUE;
1390:     TurnForward(ts);
1391:     ReCompute(ts,tjsch,e->stepnum,stepnum);
1392:     if (e->stepnum+1 == stepnum) {
1393:       StackPop(stack,&e);
1394:       ElementDestroy(stack,e);
1395:     }
1396:   } else {
1397:     /* restore a checkpoint */
1398:     StackTop(stack,&e);
1399:     UpdateTS(ts,stack,e);
1400:     /* 2 revolve actions: restore a checkpoint and then advance */
1401:     ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,tjsch->total_steps,stepnum,localstepnum,PETSC_FALSE,&store);
1402:     if (tj->monitor) {
1403:       PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);
1404:       PetscViewerASCIIPrintf(tj->monitor,"Skip the step from %D to %D (stage values already checkpointed)\n",stepnum-localstepnum+tjsch->rctx->oldcapo,stepnum-localstepnum+tjsch->rctx->oldcapo+1);
1405:       PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);
1406:     }
1407:     if (!tjsch->rctx->reverseonestep && tjsch->rctx->stepsleft > 0) tjsch->rctx->stepsleft--;
1408:     if (e->stepnum < stepnum) {
1409:       tjsch->recompute = PETSC_TRUE;
1410:       TurnForward(ts);
1411:       ReCompute(ts,tjsch,e->stepnum,stepnum);
1412:     }
1413:     if (e->stepnum == stepnum) {
1414:       StackPop(stack,&e);
1415:       ElementDestroy(stack,e);
1416:     }
1417:   }
1418:   tjsch->rctx->reverseonestep = PETSC_FALSE;
1419:   return(0);
1420: }

1422: static PetscErrorCode SetTrajRMS(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X)
1423: {
1424:   Stack          *stack = &tjsch->stack;
1425:   PetscInt       store;
1426:   StackElement   e;

1430:   if (!stack->solution_only && stepnum == 0) return(0);
1431:   if (stack->solution_only && stepnum == tjsch->total_steps) return(0);
1432:   ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,tjsch->total_steps,stepnum,stepnum,PETSC_FALSE,&store);
1433:   if (store == 1){
1434:     if (stepnum < stack->top) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_MEMC,"Illegal modification of a non-top stack element");
1435:     ElementCreate(ts,stack,&e);
1436:     ElementSet(ts,stack,&e,stepnum,time,X);
1437:     StackPush(stack,e);
1438:   } else if (store == 2) {
1439:     DumpSingle(tj,ts,stack,tjsch->rctx->check+1);
1440:   }
1441:   return(0);
1442: }

1444: static PetscErrorCode GetTrajRMS(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum)
1445: {
1446:   Stack          *stack = &tjsch->stack;
1447:   PetscInt       whattodo,shift;
1448:   PetscInt       restart;
1449:   PetscBool      ondisk;
1450:   StackElement   e;

1454:   if (stepnum == 0 || stepnum == tjsch->total_steps) {
1455:     TurnBackward(ts);
1456:     tjsch->rctx->reverseonestep = PETSC_FALSE;
1457:     return(0);
1458:   }
1459:   tjsch->rctx->capo = stepnum;
1460:   tjsch->rctx->oldcapo = tjsch->rctx->capo;
1461:   shift = 0;
1462:   whattodo = revolve_action(&tjsch->rctx->check,&tjsch->rctx->capo,&tjsch->rctx->fine,tjsch->rctx->snaps_in,&tjsch->rctx->info,&tjsch->rctx->where); /* whattodo=restore */
1463:   printwhattodo(tj->monitor,whattodo,tjsch->rctx,shift);
1464:   /* restore a checkpoint */
1465:   restart = tjsch->rctx->capo;
1466:   if (!tjsch->rctx->where) {
1467:     ondisk = PETSC_TRUE;
1468:     LoadSingle(tj,ts,stack,tjsch->rctx->check+1);
1469:     TurnBackward(ts);
1470:   } else {
1471:     ondisk = PETSC_FALSE;
1472:     StackTop(stack,&e);
1473:     UpdateTS(ts,stack,e);
1474:   }
1475:   if (!stack->solution_only) { /* whattodo must be 5 or 8 */
1476:     /* ask Revolve what to do next */
1477:     tjsch->rctx->oldcapo = tjsch->rctx->capo;
1478:     whattodo = revolve_action(&tjsch->rctx->check,&tjsch->rctx->capo,&tjsch->rctx->fine,tjsch->rctx->snaps_in,&tjsch->rctx->info,&tjsch->rctx->where); /* must return 1 or 3 or 4*/
1479:     printwhattodo(tj->monitor,whattodo,tjsch->rctx,shift);
1480:     if (whattodo == 3 || whattodo == 4) tjsch->rctx->reverseonestep = PETSC_TRUE;
1481:     if (whattodo == 1) tjsch->rctx->stepsleft = tjsch->rctx->capo-tjsch->rctx->oldcapo;
1482:     if (tj->monitor) {
1483:       PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);
1484:       PetscViewerASCIIPrintf(tj->monitor,"Skip the step from %D to %D (stage values already checkpointed)\n",tjsch->rctx->oldcapo,tjsch->rctx->oldcapo+1);
1485:       PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);
1486:     }
1487:     if (!tjsch->rctx->reverseonestep && tjsch->rctx->stepsleft > 0) tjsch->rctx->stepsleft--;
1488:     restart++; /* skip one step */
1489:   }
1490:   if (stack->solution_only || (!stack->solution_only && restart < stepnum)) {
1491:     tjsch->recompute = PETSC_TRUE;
1492:     TurnForward(ts);
1493:     ReCompute(ts,tjsch,restart,stepnum);
1494:   }
1495:   if (!ondisk && ( (stack->solution_only && e->stepnum+1 == stepnum) || (!stack->solution_only && e->stepnum == stepnum) )) {
1496:     StackPop(stack,&e);
1497:     ElementDestroy(stack,e);
1498:   }
1499:   tjsch->rctx->reverseonestep = PETSC_FALSE;
1500:   return(0);
1501: }
1502: #endif

1504: static PetscErrorCode TSTrajectorySet_Memory(TSTrajectory tj,TS ts,PetscInt stepnum,PetscReal time,Vec X)
1505: {
1506:   TJScheduler *tjsch = (TJScheduler*)tj->data;

1510:   if (!tjsch->recompute) { /* use global stepnum in the forward sweep */
1511:     TSGetStepNumber(ts,&stepnum);
1512:   }
1513:   /* for consistency */
1514:   if (!tjsch->recompute && stepnum == 0) ts->ptime_prev = ts->ptime-ts->time_step;
1515:   switch (tjsch->stype) {
1516:     case NONE:
1517:       SetTrajN(ts,tjsch,stepnum,time,X);
1518:       break;
1519:     case TWO_LEVEL_NOREVOLVE:
1520:       SetTrajTLNR(tj,ts,tjsch,stepnum,time,X);
1521:       break;
1522: #if defined(PETSC_HAVE_REVOLVE)
1523:     case TWO_LEVEL_REVOLVE:
1524:       SetTrajTLR(tj,ts,tjsch,stepnum,time,X);
1525:       break;
1526:     case TWO_LEVEL_TWO_REVOLVE:
1527:       SetTrajTLTR(tj,ts,tjsch,stepnum,time,X);
1528:       break;
1529:     case REVOLVE_OFFLINE:
1530:       SetTrajROF(tj,ts,tjsch,stepnum,time,X);
1531:       break;
1532:     case REVOLVE_ONLINE:
1533:       SetTrajRON(tj,ts,tjsch,stepnum,time,X);
1534:       break;
1535:     case REVOLVE_MULTISTAGE:
1536:       SetTrajRMS(tj,ts,tjsch,stepnum,time,X);
1537:       break;
1538: #endif
1539:     default:
1540:       break;
1541:   }
1542:   return(0);
1543: }

1545: static PetscErrorCode TSTrajectoryGet_Memory(TSTrajectory tj,TS ts,PetscInt stepnum,PetscReal *t)
1546: {
1547:   TJScheduler *tjsch = (TJScheduler*)tj->data;

1551:   if (stepnum == 0) return(0);
1552:   switch (tjsch->stype) {
1553:     case NONE:
1554:       GetTrajN(ts,tjsch,stepnum);
1555:       break;
1556:     case TWO_LEVEL_NOREVOLVE:
1557:       GetTrajTLNR(tj,ts,tjsch,stepnum);
1558:       break;
1559: #if defined(PETSC_HAVE_REVOLVE)
1560:     case TWO_LEVEL_REVOLVE:
1561:       GetTrajTLR(tj,ts,tjsch,stepnum);
1562:       break;
1563:     case TWO_LEVEL_TWO_REVOLVE:
1564:       GetTrajTLTR(tj,ts,tjsch,stepnum);
1565:       break;
1566:     case REVOLVE_OFFLINE:
1567:       GetTrajROF(tj,ts,tjsch,stepnum);
1568:       break;
1569:     case REVOLVE_ONLINE:
1570:       GetTrajRON(tj,ts,tjsch,stepnum);
1571:       break;
1572:     case REVOLVE_MULTISTAGE:
1573:       GetTrajRMS(tj,ts,tjsch,stepnum);
1574:       break;
1575: #endif
1576:     default:
1577:       break;
1578:   }
1579:   return(0);
1580: }

1582: PETSC_UNUSED static PetscErrorCode TSTrajectorySetStride_Memory(TSTrajectory tj,PetscInt stride)
1583: {
1584:   TJScheduler *tjsch = (TJScheduler*)tj->data;

1587:   tjsch->stride = stride;
1588:   return(0);
1589: }

1591: PETSC_UNUSED static PetscErrorCode TSTrajectorySetMaxCpsRAM_Memory(TSTrajectory tj,PetscInt max_cps_ram)
1592: {
1593:   TJScheduler *tjsch = (TJScheduler*)tj->data;

1596:   tjsch->max_cps_ram = max_cps_ram;
1597:   return(0);
1598: }

1600: PETSC_UNUSED static PetscErrorCode TSTrajectorySetMaxCpsDisk_Memory(TSTrajectory tj,PetscInt max_cps_disk)
1601: {
1602:   TJScheduler *tjsch = (TJScheduler*)tj->data;

1605:   tjsch->max_cps_disk = max_cps_disk;
1606:   return(0);
1607: }

1609: #if defined(PETSC_HAVE_REVOLVE)
1610: PETSC_UNUSED static PetscErrorCode TSTrajectorySetRevolveOnline(TSTrajectory tj,PetscBool use_online)
1611: {
1612:   TJScheduler *tjsch = (TJScheduler*)tj->data;

1615:   tjsch->use_online = use_online;
1616:   return(0);
1617: }
1618: #endif

1620: PETSC_UNUSED static PetscErrorCode TSTrajectorySetSaveStack(TSTrajectory tj,PetscBool save_stack)
1621: {
1622:   TJScheduler *tjsch = (TJScheduler*)tj->data;

1625:   tjsch->save_stack = save_stack;
1626:   return(0);
1627: }

1629: PETSC_UNUSED static PetscErrorCode TSTrajectorySetSolutionOnly(TSTrajectory tj,PetscBool solution_only)
1630: {
1631:   TJScheduler *tjsch = (TJScheduler*)tj->data;
1632:   Stack       *stack = &tjsch->stack;

1635:   stack->solution_only = solution_only;
1636:   return(0);
1637: }

1639: PETSC_UNUSED static PetscErrorCode TSTrajectorySetUseDRAM(TSTrajectory tj,PetscBool use_dram)
1640: {
1641:   TJScheduler *tjsch = (TJScheduler*)tj->data;

1644:   tjsch->stack.use_dram = use_dram;
1645:   return(0);
1646: }

1648: static PetscErrorCode TSTrajectorySetFromOptions_Memory(PetscOptionItems *PetscOptionsObject,TSTrajectory tj)
1649: {
1650:   TJScheduler    *tjsch = (TJScheduler*)tj->data;

1654:   PetscOptionsHead(PetscOptionsObject,"Memory based TS trajectory options");
1655:   {
1656:     PetscOptionsInt("-ts_trajectory_max_cps_ram","Maximum number of checkpoints in RAM","TSTrajectorySetMaxCpsRAM_Memory",tjsch->max_cps_ram,&tjsch->max_cps_ram,NULL);
1657:     PetscOptionsInt("-ts_trajectory_max_cps_disk","Maximum number of checkpoints on disk","TSTrajectorySetMaxCpsDisk_Memory",tjsch->max_cps_disk,&tjsch->max_cps_disk,NULL);
1658:     PetscOptionsInt("-ts_trajectory_stride","Stride to save checkpoints to file","TSTrajectorySetStride_Memory",tjsch->stride,&tjsch->stride,NULL);
1659: #if defined(PETSC_HAVE_REVOLVE)
1660:     PetscOptionsBool("-ts_trajectory_revolve_online","Trick TS trajectory into using online mode of revolve","TSTrajectorySetRevolveOnline",tjsch->use_online,&tjsch->use_online,NULL);
1661: #endif
1662:     PetscOptionsBool("-ts_trajectory_save_stack","Save all stack to disk","TSTrajectorySetSaveStack",tjsch->save_stack,&tjsch->save_stack,NULL);
1663:     PetscOptionsBool("-ts_trajectory_solution_only","Checkpoint solution only","TSTrajectorySetSolutionOnly",tjsch->stack.solution_only,&tjsch->stack.solution_only,NULL);
1664:     PetscOptionsBool("-ts_trajectory_use_dram","Use DRAM for checkpointing","TSTrajectorySetUseDRAM",tjsch->stack.use_dram,&tjsch->stack.use_dram,NULL);
1665:   }
1666:   PetscOptionsTail();
1667:   return(0);
1668: }

1670: static PetscErrorCode TSTrajectorySetUp_Memory(TSTrajectory tj,TS ts)
1671: {
1672:   TJScheduler    *tjsch = (TJScheduler*)tj->data;
1673:   Stack          *stack = &tjsch->stack;
1674: #if defined(PETSC_HAVE_REVOLVE)
1675:   RevolveCTX     *rctx,*rctx2;
1676:   DiskStack      *diskstack = &tjsch->diskstack;
1677:   PetscInt       diskblocks;
1678: #endif
1679:   PetscInt       numY;
1680:   PetscBool      flg;

1684:   PetscStrcmp(((PetscObject)ts->adapt)->type_name,TSADAPTNONE,&flg);
1685:   if (flg) tjsch->total_steps = PetscMin(ts->max_steps,(PetscInt)(PetscCeilReal((ts->max_time-ts->ptime)/ts->time_step))); /* fixed time step */
1686:   if (tjsch->max_cps_ram > 0) stack->stacksize = tjsch->max_cps_ram;

1688:   if (tjsch->stride > 1) { /* two level mode */
1689:     if (tjsch->save_stack && tjsch->max_cps_disk > 1 && tjsch->max_cps_disk <= tjsch->max_cps_ram) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_INCOMP,"The specified disk capacity is not enough to store a full stack of RAM checkpoints. You might want to change the disk capacity or use single level checkpointing instead.");
1690:     if (tjsch->max_cps_disk <= 1 && tjsch->max_cps_ram > 1 && tjsch->max_cps_ram <= tjsch->stride-1) tjsch->stype = TWO_LEVEL_REVOLVE; /* use revolve_offline for each stride */
1691:     if (tjsch->max_cps_disk > 1 && tjsch->max_cps_ram > 1 && tjsch->max_cps_ram <= tjsch->stride-1) tjsch->stype = TWO_LEVEL_TWO_REVOLVE;  /* use revolve_offline for each stride */
1692:     if (tjsch->max_cps_disk <= 1 && (tjsch->max_cps_ram >= tjsch->stride || tjsch->max_cps_ram == -1)) tjsch->stype = TWO_LEVEL_NOREVOLVE; /* can also be handled by TWO_LEVEL_REVOLVE */
1693:   } else { /* single level mode */
1694:     if (flg) { /* fixed time step */
1695:       if (tjsch->max_cps_ram >= tjsch->total_steps-1 || tjsch->max_cps_ram < 1) tjsch->stype = NONE; /* checkpoint all */
1696:       else tjsch->stype = (tjsch->max_cps_disk>1) ? REVOLVE_MULTISTAGE : REVOLVE_OFFLINE;
1697:     } else tjsch->stype = NONE; /* checkpoint all for adaptive time step */
1698: #if defined(PETSC_HAVE_REVOLVE)
1699:     if (tjsch->use_online) tjsch->stype = REVOLVE_ONLINE; /* trick into online (for testing purpose only) */
1700: #endif
1701:   }

1703:   if (tjsch->stype > TWO_LEVEL_NOREVOLVE) {
1704: #ifndef PETSC_HAVE_REVOLVE
1705:     SETERRQ(PetscObjectComm((PetscObject)ts),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.");
1706: #else
1707:     switch (tjsch->stype) {
1708:       case TWO_LEVEL_REVOLVE:
1709:         revolve_create_offline(tjsch->stride,tjsch->max_cps_ram);
1710:         break;
1711:       case TWO_LEVEL_TWO_REVOLVE:
1712:         diskblocks = tjsch->save_stack ? tjsch->max_cps_disk/(tjsch->max_cps_ram+1) : tjsch->max_cps_disk; /* The block size depends on whether the stack is saved. */
1713:         diskstack->stacksize = diskblocks;
1714:         revolve_create_offline(tjsch->stride,tjsch->max_cps_ram);
1715:         revolve2_create_offline((tjsch->total_steps+tjsch->stride-1)/tjsch->stride,diskblocks);
1716:         PetscCalloc1(1,&rctx2);
1717:         rctx2->snaps_in       = diskblocks;
1718:         rctx2->reverseonestep = PETSC_FALSE;
1719:         rctx2->check          = 0;
1720:         rctx2->oldcapo        = 0;
1721:         rctx2->capo           = 0;
1722:         rctx2->info           = 2;
1723:         rctx2->fine           = (tjsch->total_steps+tjsch->stride-1)/tjsch->stride;
1724:         tjsch->rctx2          = rctx2;
1725:         diskstack->top        = -1;
1726:         PetscMalloc1(diskstack->stacksize*sizeof(PetscInt),&diskstack->container);
1727:         break;
1728:       case REVOLVE_OFFLINE:
1729:         revolve_create_offline(tjsch->total_steps,tjsch->max_cps_ram);
1730:         break;
1731:       case REVOLVE_ONLINE:
1732:         stack->stacksize = tjsch->max_cps_ram;
1733:         revolve_create_online(tjsch->max_cps_ram);
1734:         break;
1735:       case REVOLVE_MULTISTAGE:
1736:         revolve_create_multistage(tjsch->total_steps,tjsch->max_cps_ram+tjsch->max_cps_disk,tjsch->max_cps_ram);
1737:         break;
1738:       default:
1739:         break;
1740:     }
1741:     PetscCalloc1(1,&rctx);
1742:     rctx->snaps_in       = tjsch->max_cps_ram; /* for theta methods snaps_in=2*max_cps_ram */
1743:     rctx->reverseonestep = PETSC_FALSE;
1744:     rctx->check          = 0;
1745:     rctx->oldcapo        = 0;
1746:     rctx->capo           = 0;
1747:     rctx->info           = 2;
1748:     rctx->fine           = (tjsch->stride > 1) ? tjsch->stride : tjsch->total_steps;
1749:     tjsch->rctx          = rctx;
1750:     if (tjsch->stype == REVOLVE_ONLINE) rctx->fine = -1;
1751: #endif
1752:   } else {
1753:     if (tjsch->stype == TWO_LEVEL_NOREVOLVE) stack->stacksize = tjsch->stride-1; /* need tjsch->stride-1 at most */
1754:     if (tjsch->stype == NONE) {
1755:       if (flg) stack->stacksize = stack->solution_only ? tjsch->total_steps : tjsch->total_steps-1; /* fix time step */
1756:       else { /* adaptive time step */
1757:         /* if max_cps_ram is not specified, use maximal allowed number of steps for stack size */
1758:         if(tjsch->max_cps_ram == -1) stack->stacksize = ts->max_steps < PETSC_MAX_INT ? ts->max_steps : 10000;
1759:         tjsch->total_steps = stack->solution_only ? stack->stacksize : stack->stacksize+1; /* will be updated as time integration advances */
1760:       }
1761:     }
1762:   }

1764:   tjsch->recompute = PETSC_FALSE;
1765:   TSGetStages(ts,&numY,PETSC_IGNORE);
1766:   StackCreate(stack,stack->stacksize,numY);
1767:   return(0);
1768: }

1770: static PetscErrorCode TSTrajectoryDestroy_Memory(TSTrajectory tj)
1771: {
1772:   TJScheduler    *tjsch = (TJScheduler*)tj->data;

1776:   if (tjsch->stype > TWO_LEVEL_NOREVOLVE) {
1777: #if defined(PETSC_HAVE_REVOLVE)
1778:     revolve_reset();
1779:     if (tjsch->stype == TWO_LEVEL_TWO_REVOLVE) {
1780:       revolve2_reset();
1781:       PetscFree(tjsch->diskstack.container);
1782:     }
1783: #endif
1784:   }
1785:   StackDestroy(&tjsch->stack);
1786: #if defined(PETSC_HAVE_REVOLVE)
1787:   if (tjsch->stype > TWO_LEVEL_NOREVOLVE) {
1788:     PetscFree(tjsch->rctx);
1789:     PetscFree(tjsch->rctx2);
1790:   }
1791: #endif
1792:   PetscFree(tjsch);
1793:   return(0);
1794: }

1796: /*MC
1797:       TSTRAJECTORYMEMORY - Stores each solution of the ODE/ADE in memory

1799:   Level: intermediate

1801: .seealso:  TSTrajectoryCreate(), TS, TSTrajectorySetType()

1803: M*/
1804: PETSC_EXTERN PetscErrorCode TSTrajectoryCreate_Memory(TSTrajectory tj,TS ts)
1805: {
1806:   TJScheduler    *tjsch;

1810:   tj->ops->set            = TSTrajectorySet_Memory;
1811:   tj->ops->get            = TSTrajectoryGet_Memory;
1812:   tj->ops->setup          = TSTrajectorySetUp_Memory;
1813:   tj->ops->destroy        = TSTrajectoryDestroy_Memory;
1814:   tj->ops->setfromoptions = TSTrajectorySetFromOptions_Memory;

1816:   PetscCalloc1(1,&tjsch);
1817:   tjsch->stype        = NONE;
1818:   tjsch->max_cps_ram  = -1; /* -1 indicates that it is not set */
1819:   tjsch->max_cps_disk = -1; /* -1 indicates that it is not set */
1820:   tjsch->stride       = 0; /* if not zero, two-level checkpointing will be used */
1821: #if defined(PETSC_HAVE_REVOLVE)
1822:   tjsch->use_online   = PETSC_FALSE;
1823: #endif
1824:   tjsch->save_stack   = PETSC_TRUE;

1826:   tjsch->stack.solution_only = PETSC_TRUE;

1828:   tj->data = tjsch;
1829:   return(0);
1830: }