1 #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 2 #include <petscsys.h> 3 #if defined(PETSC_HAVE_REVOLVE) 4 #include <revolve_c.h> 5 6 /* Limit Revolve to 32-bits */ 7 #define PETSC_REVOLVE_INT_MAX 2147483647 8 9 typedef int PetscRevolveInt; 10 11 PETSC_STATIC_INLINE PetscErrorCode PetscRevolveIntCast(PetscInt a,PetscRevolveInt *b) 12 { 13 PetscFunctionBegin; 14 #if defined(PETSC_USE_64BIT_INDICES) 15 *b = 0; 16 if (a > PETSC_REVOLVE_INT_MAX) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Parameter is too large for Revolve, which is restricted to 32 bit integers"); 17 #endif 18 *b = (PetscRevolveInt)(a); 19 PetscFunctionReturn(0); 20 } 21 #endif 22 #if defined(PETSC_HAVE_CAMS) 23 #include <offline_schedule.h> 24 #endif 25 26 PetscLogEvent TSTrajectory_DiskWrite, TSTrajectory_DiskRead; 27 static PetscErrorCode TSTrajectorySet_Memory(TSTrajectory,TS,PetscInt,PetscReal,Vec); 28 29 typedef enum {NONE,TWO_LEVEL_NOREVOLVE,TWO_LEVEL_REVOLVE,TWO_LEVEL_TWO_REVOLVE,REVOLVE_OFFLINE,REVOLVE_ONLINE,REVOLVE_MULTISTAGE,CAMS_OFFLINE} SchedulerType; 30 31 typedef enum {UNSET=-1,SOLUTIONONLY=0,STAGESONLY=1,SOLUTION_STAGES=2} CheckpointType; 32 33 #define HaveSolution(m) ((m) == SOLUTIONONLY || (m) == SOLUTION_STAGES) 34 #define HaveStages(m) ((m) == STAGESONLY || (m) == SOLUTION_STAGES) 35 36 typedef struct _StackElement { 37 PetscInt stepnum; 38 Vec X; 39 Vec *Y; 40 PetscReal time; 41 PetscReal timeprev; /* for no solution_only mode */ 42 PetscReal timenext; /* for solution_only mode */ 43 CheckpointType cptype; 44 } *StackElement; 45 46 #if defined(PETSC_HAVE_REVOLVE) 47 typedef struct _RevolveCTX { 48 PetscBool reverseonestep; 49 PetscRevolveInt where; 50 PetscRevolveInt snaps_in; 51 PetscRevolveInt stepsleft; 52 PetscRevolveInt check; 53 PetscRevolveInt oldcapo; 54 PetscRevolveInt capo; 55 PetscRevolveInt fine; 56 PetscRevolveInt info; 57 } RevolveCTX; 58 #endif 59 60 #if defined(PETSC_HAVE_CAMS) 61 typedef struct _CAMSCTX { 62 PetscInt lastcheckpointstep; 63 PetscInt lastcheckpointtype; 64 PetscInt num_units_avail; 65 PetscInt endstep; 66 PetscInt num_stages; 67 PetscInt nextcheckpointstep; 68 PetscInt nextcheckpointtype; /* (0) solution only (1) stages (2) solution+stages */ 69 PetscInt info; 70 } CAMSCTX; 71 #endif 72 73 typedef struct _Stack { 74 PetscInt stacksize; 75 PetscInt top; 76 StackElement *container; 77 PetscInt nallocated; 78 PetscInt numY; 79 PetscBool solution_only; 80 PetscBool use_dram; 81 } Stack; 82 83 typedef struct _DiskStack { 84 PetscInt stacksize; 85 PetscInt top; 86 PetscInt *container; 87 } DiskStack; 88 89 typedef struct _TJScheduler { 90 SchedulerType stype; 91 #if defined(PETSC_HAVE_REVOLVE) 92 RevolveCTX *rctx,*rctx2; 93 PetscBool use_online; 94 PetscInt store_stride; 95 #endif 96 #if defined(PETSC_HAVE_CAMS) 97 CAMSCTX *actx; 98 #endif 99 PetscBool recompute; 100 PetscBool skip_trajectory; 101 PetscBool save_stack; 102 PetscInt max_units_ram; /* maximum checkpointing units in RAM */ 103 PetscInt max_units_disk; /* maximum checkpointing units on disk */ 104 PetscInt max_cps_ram; /* maximum checkpoints in RAM */ 105 PetscInt max_cps_disk; /* maximum checkpoints on disk */ 106 PetscInt stride; 107 PetscInt total_steps; /* total number of steps */ 108 Stack stack; 109 DiskStack diskstack; 110 PetscViewer viewer; 111 } TJScheduler; 112 113 static PetscErrorCode TurnForwardWithStepsize(TS ts,PetscReal nextstepsize) 114 { 115 PetscErrorCode ierr; 116 117 PetscFunctionBegin; 118 /* reverse the direction */ 119 ierr = TSSetTimeStep(ts,nextstepsize);CHKERRQ(ierr); 120 PetscFunctionReturn(0); 121 } 122 123 static PetscErrorCode TurnForward(TS ts) 124 { 125 PetscReal stepsize; 126 PetscErrorCode ierr; 127 128 PetscFunctionBegin; 129 /* reverse the direction */ 130 ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr); 131 ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr); 132 PetscFunctionReturn(0); 133 } 134 135 static PetscErrorCode TurnBackward(TS ts) 136 { 137 PetscReal stepsize; 138 PetscErrorCode ierr; 139 140 PetscFunctionBegin; 141 if (!ts->trajectory->adjoint_solve_mode) PetscFunctionReturn(0); 142 /* reverse the direction */ 143 stepsize = ts->ptime_prev-ts->ptime; 144 ierr = TSSetTimeStep(ts,stepsize);CHKERRQ(ierr); 145 PetscFunctionReturn(0); 146 } 147 148 static PetscErrorCode ElementCreate(TS ts,CheckpointType cptype,Stack *stack,StackElement *e) 149 { 150 Vec X; 151 Vec *Y; 152 PetscErrorCode ierr; 153 154 PetscFunctionBegin; 155 if (stack->top < stack->stacksize-1 && stack->container[stack->top+1]) { 156 *e = stack->container[stack->top+1]; 157 if (HaveSolution(cptype) && !(*e)->X) { 158 ierr = TSGetSolution(ts,&X);CHKERRQ(ierr); 159 ierr = VecDuplicate(X,&(*e)->X);CHKERRQ(ierr); 160 } 161 if (cptype==1 && (*e)->X) { 162 ierr = VecDestroy(&(*e)->X);CHKERRQ(ierr); 163 } 164 if (HaveStages(cptype) && !(*e)->Y) { 165 ierr = TSGetStages(ts,&stack->numY,&Y);CHKERRQ(ierr); 166 if (stack->numY) { 167 ierr = VecDuplicateVecs(Y[0],stack->numY,&(*e)->Y);CHKERRQ(ierr); 168 } 169 } 170 if (cptype==0 && (*e)->Y) { 171 ierr = VecDestroyVecs(stack->numY,&(*e)->Y);CHKERRQ(ierr); 172 } 173 (*e)->cptype = cptype; 174 PetscFunctionReturn(0); 175 } 176 if (stack->use_dram) { 177 ierr = PetscMallocSetDRAM();CHKERRQ(ierr); 178 } 179 ierr = PetscNew(e);CHKERRQ(ierr); 180 if (HaveSolution(cptype)) { 181 ierr = TSGetSolution(ts,&X);CHKERRQ(ierr); 182 ierr = VecDuplicate(X,&(*e)->X);CHKERRQ(ierr); 183 } 184 if (HaveStages(cptype)) { 185 ierr = TSGetStages(ts,&stack->numY,&Y);CHKERRQ(ierr); 186 if (stack->numY) { 187 ierr = VecDuplicateVecs(Y[0],stack->numY,&(*e)->Y);CHKERRQ(ierr); 188 } 189 } 190 if (stack->use_dram) { 191 ierr = PetscMallocResetDRAM();CHKERRQ(ierr); 192 } 193 stack->nallocated++; 194 (*e)->cptype = cptype; 195 PetscFunctionReturn(0); 196 } 197 198 static PetscErrorCode ElementSet(TS ts, Stack *stack, StackElement *e, PetscInt stepnum, PetscReal time, Vec X) 199 { 200 Vec *Y; 201 PetscInt i; 202 PetscReal timeprev; 203 PetscErrorCode ierr; 204 205 PetscFunctionBegin; 206 if (HaveSolution((*e)->cptype)) { 207 ierr = VecCopy(X,(*e)->X);CHKERRQ(ierr); 208 } 209 if (HaveStages((*e)->cptype)) { 210 ierr = TSGetStages(ts,&stack->numY,&Y);CHKERRQ(ierr); 211 for (i=0;i<stack->numY;i++) { 212 ierr = VecCopy(Y[i],(*e)->Y[i]);CHKERRQ(ierr); 213 } 214 } 215 (*e)->stepnum = stepnum; 216 (*e)->time = time; 217 /* for consistency */ 218 if (stepnum == 0) { 219 (*e)->timeprev = (*e)->time - ts->time_step; 220 } else { 221 ierr = TSGetPrevTime(ts,&timeprev);CHKERRQ(ierr); 222 (*e)->timeprev = timeprev; 223 } 224 PetscFunctionReturn(0); 225 } 226 227 static PetscErrorCode ElementDestroy(Stack *stack,StackElement e) 228 { 229 PetscErrorCode ierr; 230 231 PetscFunctionBegin; 232 if (stack->use_dram) { 233 ierr = PetscMallocSetDRAM();CHKERRQ(ierr); 234 } 235 ierr = VecDestroy(&e->X);CHKERRQ(ierr); 236 if (e->Y) { 237 ierr = VecDestroyVecs(stack->numY,&e->Y);CHKERRQ(ierr); 238 } 239 ierr = PetscFree(e);CHKERRQ(ierr); 240 if (stack->use_dram) { 241 ierr = PetscMallocResetDRAM();CHKERRQ(ierr); 242 } 243 stack->nallocated--; 244 PetscFunctionReturn(0); 245 } 246 247 static PetscErrorCode StackResize(Stack *stack,PetscInt newsize) 248 { 249 StackElement *newcontainer; 250 PetscInt i; 251 PetscErrorCode ierr; 252 253 PetscFunctionBegin; 254 ierr = PetscCalloc1(newsize*sizeof(StackElement),&newcontainer);CHKERRQ(ierr); 255 for (i=0;i<stack->stacksize;i++) { 256 newcontainer[i] = stack->container[i]; 257 } 258 ierr = PetscFree(stack->container);CHKERRQ(ierr); 259 stack->container = newcontainer; 260 stack->stacksize = newsize; 261 PetscFunctionReturn(0); 262 } 263 264 static PetscErrorCode StackPush(Stack *stack,StackElement e) 265 { 266 PetscFunctionBegin; 267 if (stack->top+1 >= stack->stacksize) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MEMC,"Maximum stack size (%D) exceeded",stack->stacksize); 268 stack->container[++stack->top] = e; 269 PetscFunctionReturn(0); 270 } 271 272 static PetscErrorCode StackPop(Stack *stack,StackElement *e) 273 { 274 PetscFunctionBegin; 275 *e = NULL; 276 if (stack->top == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEMC,"Empty stack"); 277 *e = stack->container[stack->top--]; 278 PetscFunctionReturn(0); 279 } 280 281 static PetscErrorCode StackTop(Stack *stack,StackElement *e) 282 { 283 PetscFunctionBegin; 284 *e = stack->container[stack->top]; 285 PetscFunctionReturn(0); 286 } 287 288 static PetscErrorCode StackInit(Stack *stack,PetscInt size,PetscInt ny) 289 { 290 PetscErrorCode ierr; 291 292 PetscFunctionBegin; 293 stack->top = -1; 294 stack->numY = ny; 295 296 if (!stack->container) { 297 ierr = PetscCalloc1(size,&stack->container);CHKERRQ(ierr); 298 } 299 PetscFunctionReturn(0); 300 } 301 302 static PetscErrorCode StackDestroy(Stack *stack) 303 { 304 PetscInt i,n = stack->nallocated; 305 PetscErrorCode ierr; 306 307 PetscFunctionBegin; 308 if (!stack->container) PetscFunctionReturn(0); 309 if (stack->top+1 > stack->nallocated) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Stack size does not match element counter %D",stack->nallocated); 310 for (i=0; i<n; i++) { 311 ierr = ElementDestroy(stack,stack->container[i]);CHKERRQ(ierr); 312 } 313 ierr = PetscFree(stack->container);CHKERRQ(ierr); 314 PetscFunctionReturn(0); 315 } 316 317 static PetscErrorCode StackFind(Stack *stack,StackElement *e,PetscInt index) 318 { 319 PetscFunctionBegin; 320 *e = NULL; 321 if (index < 0 || index > stack->top) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Invalid index %D",index); 322 *e = stack->container[index]; 323 PetscFunctionReturn(0); 324 } 325 326 static PetscErrorCode WriteToDisk(PetscInt stepnum,PetscReal time,PetscReal timeprev,Vec X,Vec *Y,PetscInt numY,CheckpointType cptype,PetscViewer viewer) 327 { 328 PetscInt i; 329 PetscErrorCode ierr; 330 331 PetscFunctionBegin; 332 ierr = PetscViewerBinaryWrite(viewer,&stepnum,1,PETSC_INT);CHKERRQ(ierr); 333 if (HaveSolution(cptype)) { 334 ierr = VecView(X,viewer);CHKERRQ(ierr); 335 } 336 if (HaveStages(cptype)) { 337 for (i=0;i<numY;i++) { 338 ierr = VecView(Y[i],viewer);CHKERRQ(ierr); 339 } 340 } 341 ierr = PetscViewerBinaryWrite(viewer,&time,1,PETSC_REAL);CHKERRQ(ierr); 342 ierr = PetscViewerBinaryWrite(viewer,&timeprev,1,PETSC_REAL);CHKERRQ(ierr); 343 PetscFunctionReturn(0); 344 } 345 346 static PetscErrorCode ReadFromDisk(PetscInt *stepnum,PetscReal *time,PetscReal *timeprev,Vec X,Vec *Y,PetscInt numY,CheckpointType cptype,PetscViewer viewer) 347 { 348 PetscInt i; 349 PetscErrorCode ierr; 350 351 PetscFunctionBegin; 352 ierr = PetscViewerBinaryRead(viewer,stepnum,1,NULL,PETSC_INT);CHKERRQ(ierr); 353 if (HaveSolution(cptype)) { 354 ierr = VecLoad(X,viewer);CHKERRQ(ierr); 355 } 356 if (HaveStages(cptype)) { 357 for (i=0;i<numY;i++) { 358 ierr = VecLoad(Y[i],viewer);CHKERRQ(ierr); 359 } 360 } 361 ierr = PetscViewerBinaryRead(viewer,time,1,NULL,PETSC_REAL);CHKERRQ(ierr); 362 ierr = PetscViewerBinaryRead(viewer,timeprev,1,NULL,PETSC_REAL);CHKERRQ(ierr); 363 PetscFunctionReturn(0); 364 } 365 366 static PetscErrorCode StackDumpAll(TSTrajectory tj,TS ts,Stack *stack,PetscInt id) 367 { 368 Vec *Y; 369 PetscInt i,ndumped,cptype_int; 370 StackElement e = NULL; 371 TJScheduler *tjsch = (TJScheduler*)tj->data; 372 char filename[PETSC_MAX_PATH_LEN]; 373 PetscErrorCode ierr; 374 MPI_Comm comm; 375 376 PetscFunctionBegin; 377 ierr = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr); 378 if (tj->monitor) { 379 ierr = PetscViewerASCIIPushTab(tj->monitor);CHKERRQ(ierr); 380 ierr = PetscViewerASCIIPrintf(tj->monitor,"Dump stack id %D to file\n",id);CHKERRQ(ierr); 381 ierr = PetscViewerASCIIPopTab(tj->monitor);CHKERRQ(ierr); 382 } 383 ierr = PetscSNPrintf(filename,sizeof(filename),"%s/TS-STACK%06d.bin",tj->dirname,id);CHKERRQ(ierr); 384 ierr = PetscViewerFileSetName(tjsch->viewer,filename);CHKERRQ(ierr); 385 ierr = PetscViewerSetUp(tjsch->viewer);CHKERRQ(ierr); 386 ndumped = stack->top+1; 387 ierr = PetscViewerBinaryWrite(tjsch->viewer,&ndumped,1,PETSC_INT);CHKERRQ(ierr); 388 for (i=0;i<ndumped;i++) { 389 e = stack->container[i]; 390 cptype_int = (PetscInt)e->cptype; 391 ierr = PetscViewerBinaryWrite(tjsch->viewer,&cptype_int,1,PETSC_INT);CHKERRQ(ierr); 392 ierr = PetscLogEventBegin(TSTrajectory_DiskWrite,tj,ts,0,0);CHKERRQ(ierr); 393 ierr = WriteToDisk(e->stepnum,e->time,e->timeprev,e->X,e->Y,stack->numY,e->cptype,tjsch->viewer);CHKERRQ(ierr); 394 ierr = PetscLogEventEnd(TSTrajectory_DiskWrite,tj,ts,0,0);CHKERRQ(ierr); 395 ts->trajectory->diskwrites++; 396 ierr = StackPop(stack,&e);CHKERRQ(ierr); 397 } 398 /* 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 */ 399 ierr = TSGetStages(ts,&stack->numY,&Y);CHKERRQ(ierr); 400 ierr = PetscLogEventBegin(TSTrajectory_DiskWrite,tj,ts,0,0);CHKERRQ(ierr); 401 ierr = WriteToDisk(ts->steps,ts->ptime,ts->ptime_prev,ts->vec_sol,Y,stack->numY,SOLUTION_STAGES,tjsch->viewer);CHKERRQ(ierr); 402 ierr = PetscLogEventEnd(TSTrajectory_DiskWrite,tj,ts,0,0);CHKERRQ(ierr); 403 ts->trajectory->diskwrites++; 404 PetscFunctionReturn(0); 405 } 406 407 static PetscErrorCode StackLoadAll(TSTrajectory tj,TS ts,Stack *stack,PetscInt id) 408 { 409 Vec *Y; 410 PetscInt i,nloaded,cptype_int; 411 StackElement e; 412 PetscViewer viewer; 413 char filename[PETSC_MAX_PATH_LEN]; 414 PetscErrorCode ierr; 415 416 PetscFunctionBegin; 417 if (tj->monitor) { 418 ierr = PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 419 ierr = PetscViewerASCIIPrintf(tj->monitor,"Load stack from file\n");CHKERRQ(ierr); 420 ierr = PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 421 } 422 ierr = PetscSNPrintf(filename,sizeof filename,"%s/TS-STACK%06d.bin",tj->dirname,id);CHKERRQ(ierr); 423 ierr = PetscViewerBinaryOpen(PetscObjectComm((PetscObject)tj),filename,FILE_MODE_READ,&viewer);CHKERRQ(ierr); 424 ierr = PetscViewerBinarySetSkipInfo(viewer,PETSC_TRUE);CHKERRQ(ierr); 425 ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_NATIVE);CHKERRQ(ierr); 426 ierr = PetscViewerBinaryRead(viewer,&nloaded,1,NULL,PETSC_INT);CHKERRQ(ierr); 427 for (i=0;i<nloaded;i++) { 428 ierr = PetscViewerBinaryRead(viewer,&cptype_int,1,NULL,PETSC_INT);CHKERRQ(ierr); 429 ierr = ElementCreate(ts,(CheckpointType)cptype_int,stack,&e);CHKERRQ(ierr); 430 ierr = StackPush(stack,e);CHKERRQ(ierr); 431 ierr = PetscLogEventBegin(TSTrajectory_DiskRead,tj,ts,0,0);CHKERRQ(ierr); 432 ierr = ReadFromDisk(&e->stepnum,&e->time,&e->timeprev,e->X,e->Y,stack->numY,e->cptype,viewer);CHKERRQ(ierr); 433 ierr = PetscLogEventEnd(TSTrajectory_DiskRead,tj,ts,0,0);CHKERRQ(ierr); 434 ts->trajectory->diskreads++; 435 } 436 /* load the last step into TS */ 437 ierr = TSGetStages(ts,&stack->numY,&Y);CHKERRQ(ierr); 438 ierr = PetscLogEventBegin(TSTrajectory_DiskRead,tj,ts,0,0);CHKERRQ(ierr); 439 ierr = ReadFromDisk(&ts->steps,&ts->ptime,&ts->ptime_prev,ts->vec_sol,Y,stack->numY,SOLUTION_STAGES,viewer);CHKERRQ(ierr); 440 ierr = PetscLogEventEnd(TSTrajectory_DiskRead,tj,ts,0,0);CHKERRQ(ierr); 441 ts->trajectory->diskreads++; 442 ierr = TurnBackward(ts);CHKERRQ(ierr); 443 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 444 PetscFunctionReturn(0); 445 } 446 447 #if defined(PETSC_HAVE_REVOLVE) 448 static PetscErrorCode StackLoadLast(TSTrajectory tj,TS ts,Stack *stack,PetscInt id) 449 { 450 Vec *Y; 451 PetscInt size; 452 PetscViewer viewer; 453 char filename[PETSC_MAX_PATH_LEN]; 454 #if defined(PETSC_HAVE_MPIIO) 455 PetscBool usempiio; 456 #endif 457 int fd; 458 off_t off,offset; 459 PetscErrorCode ierr; 460 461 PetscFunctionBegin; 462 if (tj->monitor) { 463 ierr = PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 464 ierr = PetscViewerASCIIPrintf(tj->monitor,"Load last stack element from file\n");CHKERRQ(ierr); 465 ierr = PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 466 } 467 ierr = TSGetStages(ts,&stack->numY,&Y);CHKERRQ(ierr); 468 ierr = VecGetSize(Y[0],&size);CHKERRQ(ierr); 469 /* VecView writes to file two extra int's for class id and number of rows */ 470 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; 471 472 ierr = PetscSNPrintf(filename,sizeof filename,"%s/TS-STACK%06d.bin",tj->dirname,id);CHKERRQ(ierr); 473 ierr = PetscViewerBinaryOpen(PetscObjectComm((PetscObject)tj),filename,FILE_MODE_READ,&viewer);CHKERRQ(ierr); 474 ierr = PetscViewerBinarySetSkipInfo(viewer,PETSC_TRUE);CHKERRQ(ierr); 475 ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_NATIVE);CHKERRQ(ierr); 476 #if defined(PETSC_HAVE_MPIIO) 477 ierr = PetscViewerBinaryGetUseMPIIO(viewer,&usempiio);CHKERRQ(ierr); 478 if (usempiio) { 479 ierr = PetscViewerBinaryGetMPIIODescriptor(viewer,(MPI_File*)&fd);CHKERRQ(ierr); 480 ierr = PetscBinarySynchronizedSeek(PetscObjectComm((PetscObject)tj),fd,off,PETSC_BINARY_SEEK_END,&offset);CHKERRQ(ierr); 481 } else { 482 #endif 483 ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 484 ierr = PetscBinarySeek(fd,off,PETSC_BINARY_SEEK_END,&offset);CHKERRQ(ierr); 485 #if defined(PETSC_HAVE_MPIIO) 486 } 487 #endif 488 /* load the last step into TS */ 489 ierr = PetscLogEventBegin(TSTrajectory_DiskRead,tj,ts,0,0);CHKERRQ(ierr); 490 ierr = ReadFromDisk(&ts->steps,&ts->ptime,&ts->ptime_prev,ts->vec_sol,Y,stack->numY,SOLUTION_STAGES,viewer);CHKERRQ(ierr); 491 ierr = PetscLogEventEnd(TSTrajectory_DiskRead,tj,ts,0,0);CHKERRQ(ierr); 492 ts->trajectory->diskreads++; 493 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 494 ierr = TurnBackward(ts);CHKERRQ(ierr); 495 PetscFunctionReturn(0); 496 } 497 #endif 498 499 static PetscErrorCode DumpSingle(TSTrajectory tj,TS ts,Stack *stack,PetscInt id) 500 { 501 Vec *Y; 502 PetscInt stepnum; 503 TJScheduler *tjsch = (TJScheduler*)tj->data; 504 char filename[PETSC_MAX_PATH_LEN]; 505 PetscErrorCode ierr; 506 MPI_Comm comm; 507 508 PetscFunctionBegin; 509 ierr = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr); 510 if (tj->monitor) { 511 ierr = PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 512 ierr = PetscViewerASCIIPrintf(tj->monitor,"Dump a single point from file\n");CHKERRQ(ierr); 513 ierr = PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 514 } 515 ierr = TSGetStepNumber(ts,&stepnum);CHKERRQ(ierr); 516 ierr = PetscSNPrintf(filename,sizeof(filename),"%s/TS-CPS%06d.bin",tj->dirname,id);CHKERRQ(ierr); 517 ierr = PetscViewerFileSetName(tjsch->viewer,filename);CHKERRQ(ierr); 518 ierr = PetscViewerSetUp(tjsch->viewer);CHKERRQ(ierr); 519 520 ierr = TSGetStages(ts,&stack->numY,&Y);CHKERRQ(ierr); 521 ierr = PetscLogEventBegin(TSTrajectory_DiskWrite,tj,ts,0,0);CHKERRQ(ierr); 522 ierr = WriteToDisk(stepnum,ts->ptime,ts->ptime_prev,ts->vec_sol,Y,stack->numY,SOLUTION_STAGES,tjsch->viewer);CHKERRQ(ierr); 523 ierr = PetscLogEventEnd(TSTrajectory_DiskWrite,tj,ts,0,0);CHKERRQ(ierr); 524 ts->trajectory->diskwrites++; 525 PetscFunctionReturn(0); 526 } 527 528 static PetscErrorCode LoadSingle(TSTrajectory tj,TS ts,Stack *stack,PetscInt id) 529 { 530 Vec *Y; 531 PetscViewer viewer; 532 char filename[PETSC_MAX_PATH_LEN]; 533 PetscErrorCode ierr; 534 535 PetscFunctionBegin; 536 if (tj->monitor) { 537 ierr = PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 538 ierr = PetscViewerASCIIPrintf(tj->monitor,"Load a single point from file\n");CHKERRQ(ierr); 539 ierr = PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 540 } 541 ierr = PetscSNPrintf(filename,sizeof filename,"%s/TS-CPS%06d.bin",tj->dirname,id);CHKERRQ(ierr); 542 ierr = PetscViewerBinaryOpen(PetscObjectComm((PetscObject)tj),filename,FILE_MODE_READ,&viewer);CHKERRQ(ierr); 543 ierr = PetscViewerBinarySetSkipInfo(viewer,PETSC_TRUE);CHKERRQ(ierr); 544 ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_NATIVE);CHKERRQ(ierr); 545 ierr = TSGetStages(ts,&stack->numY,&Y);CHKERRQ(ierr); 546 ierr = PetscLogEventBegin(TSTrajectory_DiskRead,tj,ts,0,0);CHKERRQ(ierr); 547 ierr = ReadFromDisk(&ts->steps,&ts->ptime,&ts->ptime_prev,ts->vec_sol,Y,stack->numY,SOLUTION_STAGES,viewer);CHKERRQ(ierr); 548 ierr = PetscLogEventEnd(TSTrajectory_DiskRead,tj,ts,0,0);CHKERRQ(ierr); 549 ts->trajectory->diskreads++; 550 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 551 PetscFunctionReturn(0); 552 } 553 554 static PetscErrorCode UpdateTS(TS ts,Stack *stack,StackElement e,PetscInt stepnum,PetscBool adjoint_mode) 555 { 556 Vec *Y; 557 PetscInt i; 558 PetscErrorCode ierr; 559 560 PetscFunctionBegin; 561 /* In adjoint mode we do not need to copy solution if the stepnum is the same */ 562 if (!adjoint_mode || (HaveSolution(e->cptype) && e->stepnum!=stepnum)) { 563 ierr = VecCopy(e->X,ts->vec_sol);CHKERRQ(ierr); 564 } 565 if (HaveStages(e->cptype)) { 566 ierr = TSGetStages(ts,&stack->numY,&Y);CHKERRQ(ierr); 567 if (e->stepnum && e->stepnum==stepnum) { 568 for (i=0;i<stack->numY;i++) { 569 ierr = VecCopy(e->Y[i],Y[i]);CHKERRQ(ierr); 570 } 571 } else if (ts->stifflyaccurate) { 572 ierr = VecCopy(e->Y[stack->numY-1],ts->vec_sol);CHKERRQ(ierr); 573 } 574 } 575 if (adjoint_mode) { 576 ierr = TSSetTimeStep(ts,e->timeprev-e->time);CHKERRQ(ierr); /* stepsize will be negative */ 577 } else { 578 ierr = TSSetTimeStep(ts,e->time-e->timeprev);CHKERRQ(ierr); /* stepsize will be positive */ 579 } 580 ts->ptime = e->time; 581 ts->ptime_prev = e->timeprev; 582 PetscFunctionReturn(0); 583 } 584 585 static PetscErrorCode ReCompute(TS ts,TJScheduler *tjsch,PetscInt stepnumbegin,PetscInt stepnumend) 586 { 587 Stack *stack = &tjsch->stack; 588 PetscInt i; 589 PetscErrorCode ierr; 590 591 PetscFunctionBegin; 592 tjsch->recompute = PETSC_TRUE; /* hints TSTrajectorySet() that it is in recompute mode */ 593 ierr = TSSetStepNumber(ts,stepnumbegin);CHKERRQ(ierr);/* global step number */ 594 for (i=stepnumbegin;i<stepnumend;i++) { /* assume fixed step size */ 595 if (stack->solution_only && !tjsch->skip_trajectory) { /* revolve online need this */ 596 /* don't use the public interface as it will update the TSHistory: this need a better fix */ 597 ierr = TSTrajectorySet_Memory(ts->trajectory,ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 598 } 599 ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 600 ierr = TSStep(ts);CHKERRQ(ierr); 601 if (!stack->solution_only && !tjsch->skip_trajectory) { 602 /* don't use the public interface as it will update the TSHistory: this need a better fix */ 603 ierr = TSTrajectorySet_Memory(ts->trajectory,ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 604 } 605 ierr = TSEventHandler(ts);CHKERRQ(ierr); 606 if (!ts->steprollback) { 607 ierr = TSPostStep(ts);CHKERRQ(ierr); 608 } 609 } 610 ierr = TurnBackward(ts);CHKERRQ(ierr); 611 ts->trajectory->recomps += stepnumend-stepnumbegin; /* recomputation counter */ 612 ierr = TSSetStepNumber(ts,stepnumend);CHKERRQ(ierr); 613 tjsch->recompute = PETSC_FALSE; /* reset the flag for recompute mode */ 614 PetscFunctionReturn(0); 615 } 616 617 static PetscErrorCode TopLevelStore(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscInt localstepnum,PetscInt laststridesize,PetscBool *done) 618 { 619 Stack *stack = &tjsch->stack; 620 DiskStack *diskstack = &tjsch->diskstack; 621 PetscInt stridenum; 622 PetscErrorCode ierr; 623 624 PetscFunctionBegin; 625 *done = PETSC_FALSE; 626 stridenum = stepnum/tjsch->stride; 627 /* make sure saved checkpoint id starts from 1 628 skip last stride when using stridenum+1 629 skip first stride when using stridenum */ 630 if (stack->solution_only) { 631 if (tjsch->save_stack) { 632 if (localstepnum == tjsch->stride-1 && stepnum < tjsch->total_steps-laststridesize) { /* current step will be saved without going through stack */ 633 ierr = StackDumpAll(tj,ts,stack,stridenum+1);CHKERRQ(ierr); 634 if (tjsch->stype == TWO_LEVEL_TWO_REVOLVE) diskstack->container[++diskstack->top] = stridenum+1; 635 *done = PETSC_TRUE; 636 } 637 } else { 638 if (localstepnum == 0 && stepnum < tjsch->total_steps-laststridesize) { 639 ierr = DumpSingle(tj,ts,stack,stridenum+1);CHKERRQ(ierr); 640 if (tjsch->stype == TWO_LEVEL_TWO_REVOLVE) diskstack->container[++diskstack->top] = stridenum+1; 641 *done = PETSC_TRUE; 642 } 643 } 644 } else { 645 if (tjsch->save_stack) { 646 if (localstepnum == 0 && stepnum < tjsch->total_steps && stepnum != 0) { /* skip the first stride */ 647 ierr = StackDumpAll(tj,ts,stack,stridenum);CHKERRQ(ierr); 648 if (tjsch->stype == TWO_LEVEL_TWO_REVOLVE) diskstack->container[++diskstack->top] = stridenum; 649 *done = PETSC_TRUE; 650 } 651 } else { 652 if (localstepnum == 1 && stepnum < tjsch->total_steps-laststridesize) { 653 ierr = DumpSingle(tj,ts,stack,stridenum+1);CHKERRQ(ierr); 654 if (tjsch->stype == TWO_LEVEL_TWO_REVOLVE) diskstack->container[++diskstack->top] = stridenum+1; 655 *done = PETSC_TRUE; 656 } 657 } 658 } 659 PetscFunctionReturn(0); 660 } 661 662 static PetscErrorCode SetTrajN(TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X) 663 { 664 Stack *stack = &tjsch->stack; 665 StackElement e; 666 CheckpointType cptype; 667 PetscErrorCode ierr; 668 669 PetscFunctionBegin; 670 /* skip the last step */ 671 if (ts->reason) { /* only affect the forward run */ 672 /* update total_steps in the end of forward run */ 673 if (stepnum != tjsch->total_steps) tjsch->total_steps = stepnum; 674 if (stack->solution_only) { 675 /* get rid of the solution at second last step */ 676 ierr = StackPop(stack,&e);CHKERRQ(ierr); 677 } 678 PetscFunctionReturn(0); 679 } 680 /* do not save trajectory at the recompute stage for solution_only mode */ 681 if (stack->solution_only && tjsch->recompute) PetscFunctionReturn(0); 682 /* skip the first step for no_solution_only mode */ 683 if (!stack->solution_only && stepnum == 0) PetscFunctionReturn(0); 684 685 /* resize the stack */ 686 if (stack->top+1 == stack->stacksize) { 687 ierr = StackResize(stack,2*stack->stacksize);CHKERRQ(ierr); 688 } 689 /* update timenext for the previous step; necessary for step adaptivity */ 690 if (stack->top > -1) { 691 ierr = StackTop(stack,&e);CHKERRQ(ierr); 692 e->timenext = ts->ptime; 693 } 694 if (stepnum < stack->top) { 695 SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_MEMC,"Illegal modification of a non-top stack element"); 696 } 697 cptype = stack->solution_only ? SOLUTIONONLY : STAGESONLY; 698 ierr = ElementCreate(ts,cptype,stack,&e);CHKERRQ(ierr); 699 ierr = ElementSet(ts,stack,&e,stepnum,time,X);CHKERRQ(ierr); 700 ierr = StackPush(stack,e);CHKERRQ(ierr); 701 PetscFunctionReturn(0); 702 } 703 704 static PetscErrorCode SetTrajN_2(TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X) 705 { 706 Stack *stack = &tjsch->stack; 707 StackElement e; 708 CheckpointType cptype; 709 PetscErrorCode ierr; 710 711 PetscFunctionBegin; 712 if (stack->top+1 == stack->stacksize) { 713 ierr = StackResize(stack,2*stack->stacksize);CHKERRQ(ierr); 714 } 715 /* update timenext for the previous step; necessary for step adaptivity */ 716 if (stack->top > -1) { 717 ierr = StackTop(stack,&e);CHKERRQ(ierr); 718 e->timenext = ts->ptime; 719 } 720 if (stepnum < stack->top) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_MEMC,"Illegal modification of a non-top stack element"); 721 cptype = stack->solution_only ? SOLUTIONONLY : SOLUTION_STAGES; /* Always include solution in a checkpoint in non-adjoint mode */ 722 ierr = ElementCreate(ts,cptype,stack,&e);CHKERRQ(ierr); 723 ierr = ElementSet(ts,stack,&e,stepnum,time,X);CHKERRQ(ierr); 724 ierr = StackPush(stack,e);CHKERRQ(ierr); 725 PetscFunctionReturn(0); 726 } 727 728 static PetscErrorCode GetTrajN(TS ts,TJScheduler *tjsch,PetscInt stepnum) 729 { 730 Stack *stack = &tjsch->stack; 731 StackElement e; 732 PetscInt ns; 733 PetscErrorCode ierr; 734 735 PetscFunctionBegin; 736 /* If TSTrajectoryGet() is called after TSAdjointSolve() converges (e.g. outside the while loop in TSAdjointSolve()), skip getting the checkpoint. */ 737 if (ts->reason) PetscFunctionReturn(0); 738 if (stepnum == tjsch->total_steps) { 739 ierr = TurnBackward(ts);CHKERRQ(ierr); 740 PetscFunctionReturn(0); 741 } 742 /* restore a checkpoint */ 743 ierr = StackTop(stack,&e);CHKERRQ(ierr); 744 ierr = UpdateTS(ts,stack,e,stepnum,PETSC_TRUE);CHKERRQ(ierr); 745 ierr = TSGetStages(ts,&ns,PETSC_IGNORE);CHKERRQ(ierr); 746 if (stack->solution_only && ns) { /* recompute one step */ 747 ierr = TurnForwardWithStepsize(ts,e->timenext-e->time);CHKERRQ(ierr); 748 ierr = ReCompute(ts,tjsch,e->stepnum,stepnum);CHKERRQ(ierr); 749 } 750 ierr = StackPop(stack,&e);CHKERRQ(ierr); 751 PetscFunctionReturn(0); 752 } 753 754 static PetscErrorCode GetTrajN_2(TS ts,TJScheduler *tjsch,PetscInt stepnum) 755 { 756 Stack *stack = &tjsch->stack; 757 StackElement e = NULL; 758 PetscErrorCode ierr; 759 760 PetscFunctionBegin; 761 ierr = StackFind(stack,&e,stepnum);CHKERRQ(ierr); 762 if (stepnum != e->stepnum) SETERRQ2(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Inconsistent steps! %D != %D",stepnum,e->stepnum); 763 ierr = UpdateTS(ts,stack,e,stepnum,PETSC_FALSE);CHKERRQ(ierr); 764 PetscFunctionReturn(0); 765 } 766 767 static PetscErrorCode SetTrajTLNR(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X) 768 { 769 Stack *stack = &tjsch->stack; 770 PetscInt localstepnum,laststridesize; 771 StackElement e; 772 PetscBool done; 773 CheckpointType cptype; 774 PetscErrorCode ierr; 775 776 PetscFunctionBegin; 777 if (!stack->solution_only && stepnum == 0) PetscFunctionReturn(0); 778 if (stack->solution_only && stepnum == tjsch->total_steps) PetscFunctionReturn(0); 779 if (tjsch->save_stack && tjsch->recompute) PetscFunctionReturn(0); 780 781 localstepnum = stepnum%tjsch->stride; 782 /* (stridesize-1) checkpoints are saved in each stride; an extra point is added by StackDumpAll() */ 783 laststridesize = tjsch->total_steps%tjsch->stride; 784 if (!laststridesize) laststridesize = tjsch->stride; 785 786 if (!tjsch->recompute) { 787 ierr = TopLevelStore(tj,ts,tjsch,stepnum,localstepnum,laststridesize,&done);CHKERRQ(ierr); 788 if (!tjsch->save_stack && stepnum < tjsch->total_steps-laststridesize) PetscFunctionReturn(0); 789 } 790 if (!stack->solution_only && localstepnum == 0) PetscFunctionReturn(0); /* skip last point in each stride at recompute stage or last stride */ 791 if (stack->solution_only && localstepnum == tjsch->stride-1) PetscFunctionReturn(0); /* skip last step in each stride at recompute stage or last stride */ 792 793 cptype = stack->solution_only ? SOLUTIONONLY : STAGESONLY; 794 ierr = ElementCreate(ts,cptype,stack,&e);CHKERRQ(ierr); 795 ierr = ElementSet(ts,stack,&e,stepnum,time,X);CHKERRQ(ierr); 796 ierr = StackPush(stack,e);CHKERRQ(ierr); 797 PetscFunctionReturn(0); 798 } 799 800 static PetscErrorCode GetTrajTLNR(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum) 801 { 802 Stack *stack = &tjsch->stack; 803 PetscInt id,localstepnum,laststridesize; 804 StackElement e; 805 PetscErrorCode ierr; 806 807 PetscFunctionBegin; 808 if (stepnum == tjsch->total_steps) { 809 ierr = TurnBackward(ts);CHKERRQ(ierr); 810 PetscFunctionReturn(0); 811 } 812 813 localstepnum = stepnum%tjsch->stride; 814 laststridesize = tjsch->total_steps%tjsch->stride; 815 if (!laststridesize) laststridesize = tjsch->stride; 816 if (stack->solution_only) { 817 /* fill stack with info */ 818 if (localstepnum == 0 && tjsch->total_steps-stepnum >= laststridesize) { 819 id = stepnum/tjsch->stride; 820 if (tjsch->save_stack) { 821 ierr = StackLoadAll(tj,ts,stack,id);CHKERRQ(ierr); 822 tjsch->skip_trajectory = PETSC_TRUE; 823 ierr = TurnForward(ts);CHKERRQ(ierr); 824 ierr = ReCompute(ts,tjsch,id*tjsch->stride-1,id*tjsch->stride);CHKERRQ(ierr); 825 tjsch->skip_trajectory = PETSC_FALSE; 826 } else { 827 ierr = LoadSingle(tj,ts,stack,id);CHKERRQ(ierr); 828 ierr = TurnForward(ts);CHKERRQ(ierr); 829 ierr = ReCompute(ts,tjsch,(id-1)*tjsch->stride,id*tjsch->stride);CHKERRQ(ierr); 830 } 831 PetscFunctionReturn(0); 832 } 833 /* restore a checkpoint */ 834 ierr = StackPop(stack,&e);CHKERRQ(ierr); 835 ierr = UpdateTS(ts,stack,e,stepnum,PETSC_TRUE);CHKERRQ(ierr); 836 tjsch->skip_trajectory = PETSC_TRUE; 837 ierr = TurnForward(ts);CHKERRQ(ierr); 838 ierr = ReCompute(ts,tjsch,e->stepnum,stepnum);CHKERRQ(ierr); 839 tjsch->skip_trajectory = PETSC_FALSE; 840 } else { 841 CheckpointType cptype = STAGESONLY; 842 /* fill stack with info */ 843 if (localstepnum == 0 && tjsch->total_steps-stepnum >= laststridesize) { 844 id = stepnum/tjsch->stride; 845 if (tjsch->save_stack) { 846 ierr = StackLoadAll(tj,ts,stack,id);CHKERRQ(ierr); 847 } else { 848 ierr = LoadSingle(tj,ts,stack,id);CHKERRQ(ierr); 849 ierr = ElementCreate(ts,cptype,stack,&e);CHKERRQ(ierr); 850 ierr = ElementSet(ts,stack,&e,(id-1)*tjsch->stride+1,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 851 ierr = StackPush(stack,e);CHKERRQ(ierr); 852 ierr = TurnForward(ts);CHKERRQ(ierr); 853 ierr = ReCompute(ts,tjsch,e->stepnum,id*tjsch->stride);CHKERRQ(ierr); 854 } 855 PetscFunctionReturn(0); 856 } 857 /* restore a checkpoint */ 858 ierr = StackPop(stack,&e);CHKERRQ(ierr); 859 ierr = UpdateTS(ts,stack,e,stepnum,PETSC_TRUE);CHKERRQ(ierr); 860 } 861 PetscFunctionReturn(0); 862 } 863 864 #if defined(PETSC_HAVE_REVOLVE) 865 static PetscErrorCode printwhattodo(PetscViewer viewer,PetscRevolveInt whattodo,RevolveCTX *rctx,PetscRevolveInt shift) 866 { 867 PetscErrorCode ierr; 868 869 PetscFunctionBegin; 870 if (!viewer) PetscFunctionReturn(0); 871 872 switch(whattodo) { 873 case 1: 874 ierr = PetscViewerASCIIPrintf(viewer,"Advance from %D to %D\n",rctx->oldcapo+shift,rctx->capo+shift);CHKERRQ(ierr); 875 break; 876 case 2: 877 ierr = PetscViewerASCIIPrintf(viewer,"Store in checkpoint number %D (located in RAM)\n",rctx->check);CHKERRQ(ierr); 878 break; 879 case 3: 880 ierr = PetscViewerASCIIPrintf(viewer,"First turn: Initialize adjoints and reverse first step\n");CHKERRQ(ierr); 881 break; 882 case 4: 883 ierr = PetscViewerASCIIPrintf(viewer,"Forward and reverse one step\n");CHKERRQ(ierr); 884 break; 885 case 5: 886 ierr = PetscViewerASCIIPrintf(viewer,"Restore in checkpoint number %D (located in RAM)\n",rctx->check);CHKERRQ(ierr); 887 break; 888 case 7: 889 ierr = PetscViewerASCIIPrintf(viewer,"Store in checkpoint number %D (located on disk)\n",rctx->check);CHKERRQ(ierr); 890 break; 891 case 8: 892 ierr = PetscViewerASCIIPrintf(viewer,"Restore in checkpoint number %D (located on disk)\n",rctx->check);CHKERRQ(ierr); 893 break; 894 case -1: 895 ierr = PetscViewerASCIIPrintf(viewer,"Error!");CHKERRQ(ierr); 896 break; 897 } 898 PetscFunctionReturn(0); 899 } 900 901 static PetscErrorCode printwhattodo2(PetscViewer viewer,PetscRevolveInt whattodo,RevolveCTX *rctx,PetscRevolveInt shift) 902 { 903 PetscErrorCode ierr; 904 905 PetscFunctionBegin; 906 if (!viewer) PetscFunctionReturn(0); 907 908 switch(whattodo) { 909 case 1: 910 ierr = PetscViewerASCIIPrintf(viewer,"[Top Level] Advance from stride %D to stride %D\n",rctx->oldcapo+shift,rctx->capo+shift);CHKERRQ(ierr); 911 break; 912 case 2: 913 ierr = PetscViewerASCIIPrintf(viewer,"[Top Level] Store in checkpoint number %D\n",rctx->check);CHKERRQ(ierr); 914 break; 915 case 3: 916 ierr = PetscViewerASCIIPrintf(viewer,"[Top Level] First turn: Initialize adjoints and reverse first stride\n");CHKERRQ(ierr); 917 break; 918 case 4: 919 ierr = PetscViewerASCIIPrintf(viewer,"[Top Level] Forward and reverse one stride\n");CHKERRQ(ierr); 920 break; 921 case 5: 922 ierr = PetscViewerASCIIPrintf(viewer,"[Top Level] Restore in checkpoint number %D\n",rctx->check);CHKERRQ(ierr); 923 break; 924 case 7: 925 ierr = PetscViewerASCIIPrintf(viewer,"[Top Level] Store in top-level checkpoint number %D\n",rctx->check);CHKERRQ(ierr); 926 break; 927 case 8: 928 ierr = PetscViewerASCIIPrintf(viewer,"[Top Level] Restore in top-level checkpoint number %D\n",rctx->check);CHKERRQ(ierr); 929 break; 930 case -1: 931 ierr = PetscViewerASCIIPrintf(viewer,"[Top Level] Error!");CHKERRQ(ierr); 932 break; 933 } 934 PetscFunctionReturn(0); 935 } 936 937 static PetscErrorCode InitRevolve(PetscInt fine,PetscInt snaps,RevolveCTX *rctx) 938 { 939 PetscRevolveInt rsnaps,rfine; 940 PetscErrorCode ierr; 941 942 PetscFunctionBegin; 943 ierr = PetscRevolveIntCast(snaps,&rsnaps);CHKERRQ(ierr); 944 ierr = PetscRevolveIntCast(fine,&rfine);CHKERRQ(ierr); 945 revolve_reset(); 946 revolve_create_offline(rfine,rsnaps); 947 rctx->snaps_in = rsnaps; 948 rctx->fine = rfine; 949 rctx->check = 0; 950 rctx->capo = 0; 951 rctx->reverseonestep = PETSC_FALSE; 952 /* check stepsleft? */ 953 PetscFunctionReturn(0); 954 } 955 956 static PetscErrorCode FastForwardRevolve(RevolveCTX *rctx) 957 { 958 PetscRevolveInt whattodo; 959 960 PetscFunctionBegin; 961 whattodo = 0; 962 while (whattodo!=3) { /* we have to fast forward revolve to the beginning of the backward sweep due to unfriendly revolve interface */ 963 whattodo = revolve_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where); 964 } 965 PetscFunctionReturn(0); 966 } 967 968 static PetscErrorCode ApplyRevolve(PetscViewer viewer,SchedulerType stype,RevolveCTX *rctx,PetscRevolveInt total_steps,PetscRevolveInt stepnum,PetscRevolveInt localstepnum,PetscBool toplevel,PetscInt *store) 969 { 970 PetscErrorCode ierr; 971 PetscRevolveInt shift,whattodo; 972 973 PetscFunctionBegin; 974 *store = 0; 975 if (rctx->stepsleft > 0) { /* advance the solution without checkpointing anything as Revolve requires */ 976 rctx->stepsleft--; 977 PetscFunctionReturn(0); 978 } 979 /* let Revolve determine what to do next */ 980 shift = stepnum-localstepnum; 981 rctx->oldcapo = rctx->capo; 982 rctx->capo = localstepnum; 983 984 if (!toplevel) whattodo = revolve_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where); 985 else whattodo = revolve2_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where); 986 if (stype == REVOLVE_ONLINE && whattodo == 8) whattodo = 5; 987 if (stype == REVOLVE_ONLINE && whattodo == 7) whattodo = 2; 988 if (!toplevel) {ierr = printwhattodo(viewer,whattodo,rctx,shift);CHKERRQ(ierr);} 989 else {ierr = printwhattodo2(viewer,whattodo,rctx,shift);CHKERRQ(ierr);} 990 if (whattodo == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in the Revolve library"); 991 if (whattodo == 1) { /* advance some time steps */ 992 if (stype == REVOLVE_ONLINE && rctx->capo >= total_steps-1) { 993 revolve_turn(total_steps,&rctx->capo,&rctx->fine); 994 if (!toplevel) {ierr = printwhattodo(viewer,whattodo,rctx,shift);CHKERRQ(ierr);} 995 else {ierr = printwhattodo2(viewer,whattodo,rctx,shift);CHKERRQ(ierr);} 996 } 997 rctx->stepsleft = rctx->capo-rctx->oldcapo-1; 998 } 999 if (whattodo == 3 || whattodo == 4) { /* ready for a reverse step */ 1000 rctx->reverseonestep = PETSC_TRUE; 1001 } 1002 if (whattodo == 5) { /* restore a checkpoint and ask Revolve what to do next */ 1003 rctx->oldcapo = rctx->capo; 1004 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*/ 1005 else whattodo = revolve2_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where); 1006 if (!toplevel) {ierr = printwhattodo(viewer,whattodo,rctx,shift);CHKERRQ(ierr);} 1007 else {ierr = printwhattodo2(viewer,whattodo,rctx,shift);CHKERRQ(ierr);} 1008 if (whattodo == 3 || whattodo == 4) rctx->reverseonestep = PETSC_TRUE; 1009 if (whattodo == 1) rctx->stepsleft = rctx->capo-rctx->oldcapo; 1010 } 1011 if (whattodo == 7) { /* save the checkpoint to disk */ 1012 *store = 2; 1013 rctx->oldcapo = rctx->capo; 1014 whattodo = revolve_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where); /* must return 1 */ 1015 ierr = printwhattodo(viewer,whattodo,rctx,shift);CHKERRQ(ierr); 1016 rctx->stepsleft = rctx->capo-rctx->oldcapo-1; 1017 } 1018 if (whattodo == 2) { /* store a checkpoint to RAM and ask Revolve how many time steps to advance next */ 1019 *store = 1; 1020 rctx->oldcapo = rctx->capo; 1021 if (!toplevel) whattodo = revolve_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where); /* must return 1 */ 1022 else whattodo = revolve2_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where); 1023 if (!toplevel) {ierr = printwhattodo(viewer,whattodo,rctx,shift);CHKERRQ(ierr);} 1024 else {ierr = printwhattodo2(viewer,whattodo,rctx,shift);CHKERRQ(ierr);} 1025 if (stype == REVOLVE_ONLINE && rctx->capo >= total_steps-1) { 1026 revolve_turn(total_steps,&rctx->capo,&rctx->fine); 1027 ierr = printwhattodo(viewer,whattodo,rctx,shift);CHKERRQ(ierr); 1028 } 1029 rctx->stepsleft = rctx->capo-rctx->oldcapo-1; 1030 } 1031 PetscFunctionReturn(0); 1032 } 1033 1034 static PetscErrorCode SetTrajROF(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X) 1035 { 1036 Stack *stack = &tjsch->stack; 1037 PetscInt store; 1038 StackElement e; 1039 PetscRevolveInt rtotal_steps,rstepnum; 1040 CheckpointType cptype; 1041 PetscErrorCode ierr; 1042 1043 PetscFunctionBegin; 1044 if (!stack->solution_only && stepnum == 0) PetscFunctionReturn(0); 1045 if (stack->solution_only && stepnum == tjsch->total_steps) PetscFunctionReturn(0); 1046 ierr = PetscRevolveIntCast(tjsch->total_steps,&rtotal_steps);CHKERRQ(ierr); 1047 ierr = PetscRevolveIntCast(stepnum,&rstepnum);CHKERRQ(ierr); 1048 ierr = ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,rtotal_steps,rstepnum,rstepnum,PETSC_FALSE,&store);CHKERRQ(ierr); 1049 if (store == 1) { 1050 if (stepnum < stack->top) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_MEMC,"Illegal modification of a non-top stack element"); 1051 cptype = stack->solution_only ? SOLUTIONONLY : SOLUTION_STAGES; 1052 ierr = ElementCreate(ts,cptype,stack,&e);CHKERRQ(ierr); 1053 ierr = ElementSet(ts,stack,&e,stepnum,time,X);CHKERRQ(ierr); 1054 ierr = StackPush(stack,e);CHKERRQ(ierr); 1055 } 1056 PetscFunctionReturn(0); 1057 } 1058 1059 static PetscErrorCode GetTrajROF(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum) 1060 { 1061 Stack *stack = &tjsch->stack; 1062 PetscInt store; 1063 PetscRevolveInt whattodo,shift,rtotal_steps,rstepnum; 1064 StackElement e; 1065 PetscErrorCode ierr; 1066 1067 PetscFunctionBegin; 1068 if (stepnum == 0 || stepnum == tjsch->total_steps) { 1069 ierr = TurnBackward(ts);CHKERRQ(ierr); 1070 tjsch->rctx->reverseonestep = PETSC_FALSE; 1071 PetscFunctionReturn(0); 1072 } 1073 /* restore a checkpoint */ 1074 ierr = StackTop(stack,&e);CHKERRQ(ierr); 1075 ierr = UpdateTS(ts,stack,e,stepnum,PETSC_TRUE);CHKERRQ(ierr); 1076 ierr = PetscRevolveIntCast(tjsch->total_steps,&rtotal_steps);CHKERRQ(ierr); 1077 ierr = PetscRevolveIntCast(stepnum,&rstepnum);CHKERRQ(ierr); 1078 if (stack->solution_only) { /* start with restoring a checkpoint */ 1079 tjsch->rctx->capo = rstepnum; 1080 tjsch->rctx->oldcapo = tjsch->rctx->capo; 1081 shift = 0; 1082 whattodo = revolve_action(&tjsch->rctx->check,&tjsch->rctx->capo,&tjsch->rctx->fine,tjsch->rctx->snaps_in,&tjsch->rctx->info,&tjsch->rctx->where); 1083 ierr = printwhattodo(tj->monitor,whattodo,tjsch->rctx,shift);CHKERRQ(ierr); 1084 } else { /* 2 revolve actions: restore a checkpoint and then advance */ 1085 ierr = ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,rtotal_steps,rstepnum,rstepnum,PETSC_FALSE,&store);CHKERRQ(ierr); 1086 if (tj->monitor) { 1087 ierr = PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 1088 ierr = PetscViewerASCIIPrintf(tj->monitor,"Skip the step from %D to %D (stage values already checkpointed)\n",tjsch->rctx->oldcapo,tjsch->rctx->oldcapo+1);CHKERRQ(ierr); 1089 ierr = PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 1090 } 1091 if (!tjsch->rctx->reverseonestep && tjsch->rctx->stepsleft > 0) tjsch->rctx->stepsleft--; 1092 } 1093 if (stack->solution_only || (!stack->solution_only && e->stepnum < stepnum)) { 1094 ierr = TurnForward(ts);CHKERRQ(ierr); 1095 ierr = ReCompute(ts,tjsch,e->stepnum,stepnum);CHKERRQ(ierr); 1096 } 1097 if ((stack->solution_only && e->stepnum+1 == stepnum) || (!stack->solution_only && e->stepnum == stepnum)) { 1098 ierr = StackPop(stack,&e);CHKERRQ(ierr); 1099 } 1100 tjsch->rctx->reverseonestep = PETSC_FALSE; 1101 PetscFunctionReturn(0); 1102 } 1103 1104 static PetscErrorCode SetTrajRON(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X) 1105 { 1106 Stack *stack = &tjsch->stack; 1107 Vec *Y; 1108 PetscInt i,store; 1109 PetscReal timeprev; 1110 StackElement e; 1111 RevolveCTX *rctx = tjsch->rctx; 1112 PetscRevolveInt rtotal_steps,rstepnum; 1113 CheckpointType cptype; 1114 PetscErrorCode ierr; 1115 1116 PetscFunctionBegin; 1117 if (!stack->solution_only && stepnum == 0) PetscFunctionReturn(0); 1118 if (stack->solution_only && stepnum == tjsch->total_steps) PetscFunctionReturn(0); 1119 ierr = PetscRevolveIntCast(tjsch->total_steps,&rtotal_steps);CHKERRQ(ierr); 1120 ierr = PetscRevolveIntCast(stepnum,&rstepnum);CHKERRQ(ierr); 1121 ierr = ApplyRevolve(tj->monitor,tjsch->stype,rctx,rtotal_steps,rstepnum,rstepnum,PETSC_FALSE,&store);CHKERRQ(ierr); 1122 if (store == 1) { 1123 if (rctx->check != stack->top+1) { /* overwrite some non-top checkpoint in the stack */ 1124 ierr = StackFind(stack,&e,rctx->check);CHKERRQ(ierr); 1125 if (HaveSolution(e->cptype)) { 1126 ierr = VecCopy(X,e->X);CHKERRQ(ierr); 1127 } 1128 if (HaveStages(e->cptype)) { 1129 ierr = TSGetStages(ts,&stack->numY,&Y);CHKERRQ(ierr); 1130 for (i=0;i<stack->numY;i++) { 1131 ierr = VecCopy(Y[i],e->Y[i]);CHKERRQ(ierr); 1132 } 1133 } 1134 e->stepnum = stepnum; 1135 e->time = time; 1136 ierr = TSGetPrevTime(ts,&timeprev);CHKERRQ(ierr); 1137 e->timeprev = timeprev; 1138 } else { 1139 if (stepnum < stack->top) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_MEMC,"Illegal modification of a non-top stack element"); 1140 cptype = stack->solution_only ? SOLUTIONONLY : SOLUTION_STAGES; 1141 ierr = ElementCreate(ts,cptype,stack,&e);CHKERRQ(ierr); 1142 ierr = ElementSet(ts,stack,&e,stepnum,time,X);CHKERRQ(ierr); 1143 ierr = StackPush(stack,e);CHKERRQ(ierr); 1144 } 1145 } 1146 PetscFunctionReturn(0); 1147 } 1148 1149 static PetscErrorCode GetTrajRON(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum) 1150 { 1151 Stack *stack = &tjsch->stack; 1152 PetscRevolveInt whattodo,shift,rstepnum; 1153 StackElement e; 1154 PetscErrorCode ierr; 1155 1156 PetscFunctionBegin; 1157 if (stepnum == 0 || stepnum == tjsch->total_steps) { 1158 ierr = TurnBackward(ts);CHKERRQ(ierr); 1159 tjsch->rctx->reverseonestep = PETSC_FALSE; 1160 PetscFunctionReturn(0); 1161 } 1162 ierr = PetscRevolveIntCast(stepnum,&rstepnum);CHKERRQ(ierr); 1163 tjsch->rctx->capo = rstepnum; 1164 tjsch->rctx->oldcapo = tjsch->rctx->capo; 1165 shift = 0; 1166 whattodo = revolve_action(&tjsch->rctx->check,&tjsch->rctx->capo,&tjsch->rctx->fine,tjsch->rctx->snaps_in,&tjsch->rctx->info,&tjsch->rctx->where); /* whattodo=restore */ 1167 if (whattodo == 8) whattodo = 5; 1168 ierr = printwhattodo(tj->monitor,whattodo,tjsch->rctx,shift);CHKERRQ(ierr); 1169 /* restore a checkpoint */ 1170 ierr = StackFind(stack,&e,tjsch->rctx->check);CHKERRQ(ierr); 1171 ierr = UpdateTS(ts,stack,e,stepnum,PETSC_TRUE);CHKERRQ(ierr); 1172 if (!stack->solution_only) { /* whattodo must be 5 */ 1173 /* ask Revolve what to do next */ 1174 tjsch->rctx->oldcapo = tjsch->rctx->capo; 1175 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*/ 1176 ierr = printwhattodo(tj->monitor,whattodo,tjsch->rctx,shift);CHKERRQ(ierr); 1177 if (whattodo == 3 || whattodo == 4) tjsch->rctx->reverseonestep = PETSC_TRUE; 1178 if (whattodo == 1) tjsch->rctx->stepsleft = tjsch->rctx->capo-tjsch->rctx->oldcapo; 1179 if (tj->monitor) { 1180 ierr = PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 1181 ierr = PetscViewerASCIIPrintf(tj->monitor,"Skip the step from %D to %D (stage values already checkpointed)\n",tjsch->rctx->oldcapo,tjsch->rctx->oldcapo+1);CHKERRQ(ierr); 1182 ierr = PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 1183 } 1184 if (!tjsch->rctx->reverseonestep && tjsch->rctx->stepsleft > 0) tjsch->rctx->stepsleft--; 1185 } 1186 if (stack->solution_only || (!stack->solution_only && e->stepnum < stepnum)) { 1187 ierr = TurnForward(ts);CHKERRQ(ierr); 1188 ierr = ReCompute(ts,tjsch,e->stepnum,stepnum);CHKERRQ(ierr); 1189 } 1190 tjsch->rctx->reverseonestep = PETSC_FALSE; 1191 PetscFunctionReturn(0); 1192 } 1193 1194 static PetscErrorCode SetTrajTLR(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X) 1195 { 1196 Stack *stack = &tjsch->stack; 1197 PetscInt store,localstepnum,laststridesize; 1198 StackElement e; 1199 PetscBool done = PETSC_FALSE; 1200 PetscRevolveInt rtotal_steps,rstepnum,rlocalstepnum; 1201 CheckpointType cptype; 1202 PetscErrorCode ierr; 1203 1204 PetscFunctionBegin; 1205 if (!stack->solution_only && stepnum == 0) PetscFunctionReturn(0); 1206 if (stack->solution_only && stepnum == tjsch->total_steps) PetscFunctionReturn(0); 1207 1208 localstepnum = stepnum%tjsch->stride; 1209 laststridesize = tjsch->total_steps%tjsch->stride; 1210 if (!laststridesize) laststridesize = tjsch->stride; 1211 1212 if (!tjsch->recompute) { 1213 ierr = TopLevelStore(tj,ts,tjsch,stepnum,localstepnum,laststridesize,&done);CHKERRQ(ierr); 1214 /* revolve is needed for the last stride; different starting points for last stride between solutin_only and !solutin_only */ 1215 if (!stack->solution_only && !tjsch->save_stack && stepnum <= tjsch->total_steps-laststridesize) PetscFunctionReturn(0); 1216 if (stack->solution_only && !tjsch->save_stack && stepnum < tjsch->total_steps-laststridesize) PetscFunctionReturn(0); 1217 } 1218 if (tjsch->save_stack && done) { 1219 ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr); 1220 PetscFunctionReturn(0); 1221 } 1222 if (laststridesize < tjsch->stride) { 1223 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 */ 1224 ierr = InitRevolve(laststridesize,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr); 1225 } 1226 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 */ 1227 ierr = InitRevolve(laststridesize,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr); 1228 } 1229 } 1230 ierr = PetscRevolveIntCast(tjsch->total_steps,&rtotal_steps);CHKERRQ(ierr); 1231 ierr = PetscRevolveIntCast(stepnum,&rstepnum);CHKERRQ(ierr); 1232 ierr = PetscRevolveIntCast(localstepnum,&rlocalstepnum);CHKERRQ(ierr); 1233 ierr = ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,rtotal_steps,rstepnum,rlocalstepnum,PETSC_FALSE,&store);CHKERRQ(ierr); 1234 if (store == 1) { 1235 if (localstepnum < stack->top) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_MEMC,"Illegal modification of a non-top stack element"); 1236 cptype = stack->solution_only ? SOLUTIONONLY : SOLUTION_STAGES; 1237 ierr = ElementCreate(ts,cptype,stack,&e);CHKERRQ(ierr); 1238 ierr = ElementSet(ts,stack,&e,stepnum,time,X);CHKERRQ(ierr); 1239 ierr = StackPush(stack,e);CHKERRQ(ierr); 1240 } 1241 PetscFunctionReturn(0); 1242 } 1243 1244 static PetscErrorCode GetTrajTLR(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum) 1245 { 1246 Stack *stack = &tjsch->stack; 1247 PetscRevolveInt whattodo,shift,rstepnum,rlocalstepnum,rtotal_steps; 1248 PetscInt localstepnum,stridenum,laststridesize,store; 1249 StackElement e; 1250 CheckpointType cptype; 1251 PetscErrorCode ierr; 1252 1253 PetscFunctionBegin; 1254 localstepnum = stepnum%tjsch->stride; 1255 stridenum = stepnum/tjsch->stride; 1256 if (stepnum == tjsch->total_steps) { 1257 ierr = TurnBackward(ts);CHKERRQ(ierr); 1258 tjsch->rctx->reverseonestep = PETSC_FALSE; 1259 PetscFunctionReturn(0); 1260 } 1261 laststridesize = tjsch->total_steps%tjsch->stride; 1262 if (!laststridesize) laststridesize = tjsch->stride; 1263 ierr = PetscRevolveIntCast(tjsch->total_steps,&rtotal_steps);CHKERRQ(ierr); 1264 ierr = PetscRevolveIntCast(stepnum,&rstepnum);CHKERRQ(ierr); 1265 ierr = PetscRevolveIntCast(localstepnum,&rlocalstepnum);CHKERRQ(ierr); 1266 if (stack->solution_only) { 1267 /* fill stack */ 1268 if (localstepnum == 0 && stepnum <= tjsch->total_steps-laststridesize) { 1269 if (tjsch->save_stack) { 1270 ierr = StackLoadAll(tj,ts,stack,stridenum);CHKERRQ(ierr); 1271 ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr); 1272 ierr = FastForwardRevolve(tjsch->rctx);CHKERRQ(ierr); 1273 tjsch->skip_trajectory = PETSC_TRUE; 1274 ierr = TurnForward(ts);CHKERRQ(ierr); 1275 ierr = ReCompute(ts,tjsch,stridenum*tjsch->stride-1,stridenum*tjsch->stride);CHKERRQ(ierr); 1276 tjsch->skip_trajectory = PETSC_FALSE; 1277 } else { 1278 ierr = LoadSingle(tj,ts,stack,stridenum);CHKERRQ(ierr); 1279 ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr); 1280 ierr = TurnForward(ts);CHKERRQ(ierr); 1281 ierr = ReCompute(ts,tjsch,(stridenum-1)*tjsch->stride,stridenum*tjsch->stride);CHKERRQ(ierr); 1282 } 1283 PetscFunctionReturn(0); 1284 } 1285 /* restore a checkpoint */ 1286 ierr = StackTop(stack,&e);CHKERRQ(ierr); 1287 ierr = UpdateTS(ts,stack,e,stepnum,PETSC_TRUE);CHKERRQ(ierr); 1288 /* start with restoring a checkpoint */ 1289 tjsch->rctx->capo = rstepnum; 1290 tjsch->rctx->oldcapo = tjsch->rctx->capo; 1291 shift = rstepnum-rlocalstepnum; 1292 whattodo = revolve_action(&tjsch->rctx->check,&tjsch->rctx->capo,&tjsch->rctx->fine,tjsch->rctx->snaps_in,&tjsch->rctx->info,&tjsch->rctx->where); 1293 ierr = printwhattodo(tj->monitor,whattodo,tjsch->rctx,shift);CHKERRQ(ierr); 1294 ierr = TurnForward(ts);CHKERRQ(ierr); 1295 ierr = ReCompute(ts,tjsch,e->stepnum,stepnum);CHKERRQ(ierr); 1296 if (e->stepnum+1 == stepnum) { 1297 ierr = StackPop(stack,&e);CHKERRQ(ierr); 1298 } 1299 } else { 1300 /* fill stack with info */ 1301 if (localstepnum == 0 && tjsch->total_steps-stepnum >= laststridesize) { 1302 if (tjsch->save_stack) { 1303 ierr = StackLoadAll(tj,ts,stack,stridenum);CHKERRQ(ierr); 1304 ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr); 1305 ierr = FastForwardRevolve(tjsch->rctx);CHKERRQ(ierr); 1306 } else { 1307 PetscRevolveInt rnum; 1308 ierr = LoadSingle(tj,ts,stack,stridenum);CHKERRQ(ierr); 1309 ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr); 1310 ierr = PetscRevolveIntCast((stridenum-1)*tjsch->stride+1,&rnum);CHKERRQ(ierr); 1311 ierr = ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,rtotal_steps,rnum,1,PETSC_FALSE,&store);CHKERRQ(ierr); 1312 if (tj->monitor) { 1313 ierr = PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 1314 ierr = 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);CHKERRQ(ierr); 1315 ierr = PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 1316 } 1317 cptype = SOLUTION_STAGES; 1318 ierr = ElementCreate(ts,cptype,stack,&e);CHKERRQ(ierr); 1319 ierr = ElementSet(ts,stack,&e,(stridenum-1)*tjsch->stride+1,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 1320 ierr = StackPush(stack,e);CHKERRQ(ierr); 1321 ierr = TurnForward(ts);CHKERRQ(ierr); 1322 ierr = ReCompute(ts,tjsch,e->stepnum,stridenum*tjsch->stride);CHKERRQ(ierr); 1323 } 1324 PetscFunctionReturn(0); 1325 } 1326 /* restore a checkpoint */ 1327 ierr = StackTop(stack,&e);CHKERRQ(ierr); 1328 ierr = UpdateTS(ts,stack,e,stepnum,PETSC_TRUE);CHKERRQ(ierr); 1329 /* 2 revolve actions: restore a checkpoint and then advance */ 1330 ierr = ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,rtotal_steps,rstepnum,rlocalstepnum,PETSC_FALSE,&store);CHKERRQ(ierr); 1331 if (tj->monitor) { 1332 ierr = PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 1333 ierr = 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);CHKERRQ(ierr); 1334 ierr = PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 1335 } 1336 if (!tjsch->rctx->reverseonestep && tjsch->rctx->stepsleft > 0) tjsch->rctx->stepsleft--; 1337 if (e->stepnum < stepnum) { 1338 ierr = TurnForward(ts);CHKERRQ(ierr); 1339 ierr = ReCompute(ts,tjsch,e->stepnum,stepnum);CHKERRQ(ierr); 1340 } 1341 if (e->stepnum == stepnum) { 1342 ierr = StackPop(stack,&e);CHKERRQ(ierr); 1343 } 1344 } 1345 tjsch->rctx->reverseonestep = PETSC_FALSE; 1346 PetscFunctionReturn(0); 1347 } 1348 1349 static PetscErrorCode SetTrajTLTR(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X) 1350 { 1351 Stack *stack = &tjsch->stack; 1352 PetscInt store,localstepnum,stridenum,laststridesize; 1353 StackElement e; 1354 PetscBool done = PETSC_FALSE; 1355 PetscRevolveInt rlocalstepnum,rstepnum,rtotal_steps; 1356 PetscErrorCode ierr; 1357 1358 PetscFunctionBegin; 1359 if (!stack->solution_only && stepnum == 0) PetscFunctionReturn(0); 1360 if (stack->solution_only && stepnum == tjsch->total_steps) PetscFunctionReturn(0); 1361 1362 localstepnum = stepnum%tjsch->stride; /* index at the bottom level (inside a stride) */ 1363 stridenum = stepnum/tjsch->stride; /* index at the top level */ 1364 laststridesize = tjsch->total_steps%tjsch->stride; 1365 if (!laststridesize) laststridesize = tjsch->stride; 1366 if (stack->solution_only && localstepnum == 0 && !tjsch->rctx2->reverseonestep) { 1367 ierr = PetscRevolveIntCast((tjsch->total_steps+tjsch->stride-1)/tjsch->stride,&rtotal_steps);CHKERRQ(ierr); 1368 ierr = PetscRevolveIntCast(stridenum,&rstepnum);CHKERRQ(ierr); 1369 ierr = ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx2,rtotal_steps,rstepnum,rstepnum,PETSC_TRUE,&tjsch->store_stride);CHKERRQ(ierr); 1370 if (laststridesize < tjsch->stride && stepnum == tjsch->total_steps-laststridesize) { 1371 ierr = InitRevolve(laststridesize,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr); 1372 } 1373 } 1374 if (!stack->solution_only && localstepnum == 1 && !tjsch->rctx2->reverseonestep) { 1375 ierr = PetscRevolveIntCast((tjsch->total_steps+tjsch->stride-1)/tjsch->stride,&rtotal_steps);CHKERRQ(ierr); 1376 ierr = PetscRevolveIntCast(stridenum,&rstepnum);CHKERRQ(ierr); 1377 ierr = ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx2,rtotal_steps,rstepnum,rstepnum,PETSC_TRUE,&tjsch->store_stride);CHKERRQ(ierr); 1378 if (laststridesize < tjsch->stride && stepnum == tjsch->total_steps-laststridesize+1) { 1379 ierr = InitRevolve(laststridesize,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr); 1380 } 1381 } 1382 if (tjsch->store_stride) { 1383 ierr = TopLevelStore(tj,ts,tjsch,stepnum,localstepnum,laststridesize,&done);CHKERRQ(ierr); 1384 if (done) { 1385 ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr); 1386 PetscFunctionReturn(0); 1387 } 1388 } 1389 if (stepnum < tjsch->total_steps-laststridesize) { 1390 if (tjsch->save_stack && !tjsch->store_stride && !tjsch->rctx2->reverseonestep) PetscFunctionReturn(0); /* store or forward-and-reverse at top level trigger revolve at bottom level */ 1391 if (!tjsch->save_stack && !tjsch->rctx2->reverseonestep) PetscFunctionReturn(0); /* store operation does not require revolve be called at bottom level */ 1392 } 1393 /* 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 */ 1394 if (!stack->solution_only && localstepnum == 0 && stepnum != tjsch->total_steps && !tjsch->recompute) PetscFunctionReturn(0); 1395 ierr = PetscRevolveIntCast(tjsch->total_steps,&rtotal_steps);CHKERRQ(ierr); 1396 ierr = PetscRevolveIntCast(stepnum,&rstepnum);CHKERRQ(ierr); 1397 ierr = PetscRevolveIntCast(localstepnum,&rlocalstepnum);CHKERRQ(ierr); 1398 ierr = ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,rtotal_steps,rstepnum,rlocalstepnum,PETSC_FALSE,&store);CHKERRQ(ierr); 1399 if (store == 1) { 1400 CheckpointType cptype; 1401 if (localstepnum < stack->top) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_MEMC,"Illegal modification of a non-top stack element"); 1402 cptype = stack->solution_only ? SOLUTIONONLY : SOLUTION_STAGES; 1403 ierr = ElementCreate(ts,cptype,stack,&e);CHKERRQ(ierr); 1404 ierr = ElementSet(ts,stack,&e,stepnum,time,X);CHKERRQ(ierr); 1405 ierr = StackPush(stack,e);CHKERRQ(ierr); 1406 } 1407 PetscFunctionReturn(0); 1408 } 1409 1410 static PetscErrorCode GetTrajTLTR(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum) 1411 { 1412 Stack *stack = &tjsch->stack; 1413 DiskStack *diskstack = &tjsch->diskstack; 1414 PetscInt localstepnum,stridenum,restoredstridenum,laststridesize,store; 1415 StackElement e; 1416 PetscRevolveInt whattodo,shift; 1417 PetscRevolveInt rtotal_steps,rstepnum,rlocalstepnum; 1418 PetscErrorCode ierr; 1419 1420 PetscFunctionBegin; 1421 localstepnum = stepnum%tjsch->stride; 1422 stridenum = stepnum/tjsch->stride; 1423 if (stepnum == tjsch->total_steps) { 1424 ierr = TurnBackward(ts);CHKERRQ(ierr); 1425 tjsch->rctx->reverseonestep = PETSC_FALSE; 1426 PetscFunctionReturn(0); 1427 } 1428 laststridesize = tjsch->total_steps%tjsch->stride; 1429 if (!laststridesize) laststridesize = tjsch->stride; 1430 /* 1431 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: 1432 Case 1 (save_stack) 1433 Restore a disk checkpoint; update TS with the last element in the restored data; recompute to the current point. 1434 Case 2 (!save_stack) 1435 Restore a disk checkpoint; update TS with the restored point; recompute to the current point. 1436 */ 1437 if (localstepnum == 0 && stepnum <= tjsch->total_steps-laststridesize) { 1438 /* restore the top element in the stack for disk checkpoints */ 1439 restoredstridenum = diskstack->container[diskstack->top]; 1440 tjsch->rctx2->reverseonestep = PETSC_FALSE; 1441 /* top-level revolve must be applied before current step, just like the solution_only mode for single-level revolve */ 1442 if (!tjsch->save_stack && stack->solution_only) { /* start with restoring a checkpoint */ 1443 ierr = PetscRevolveIntCast(stridenum,&rstepnum);CHKERRQ(ierr); 1444 tjsch->rctx2->capo = rstepnum; 1445 tjsch->rctx2->oldcapo = tjsch->rctx2->capo; 1446 shift = 0; 1447 whattodo = revolve2_action(&tjsch->rctx2->check,&tjsch->rctx2->capo,&tjsch->rctx2->fine,tjsch->rctx2->snaps_in,&tjsch->rctx2->info,&tjsch->rctx2->where); 1448 ierr = printwhattodo2(tj->monitor,whattodo,tjsch->rctx2,shift);CHKERRQ(ierr); 1449 } else { /* 2 revolve actions: restore a checkpoint and then advance */ 1450 ierr = PetscRevolveIntCast((tjsch->total_steps+tjsch->stride-1)/tjsch->stride,&rtotal_steps);CHKERRQ(ierr); 1451 ierr = PetscRevolveIntCast(stridenum,&rstepnum);CHKERRQ(ierr); 1452 ierr = ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx2,rtotal_steps,rstepnum,rstepnum,PETSC_TRUE,&tjsch->store_stride);CHKERRQ(ierr); 1453 if (tj->monitor) { 1454 ierr = PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 1455 ierr = PetscViewerASCIIPrintf(tj->monitor,"[Top Level] Skip the stride from %D to %D (stage values already checkpointed)\n",tjsch->rctx2->oldcapo,tjsch->rctx2->oldcapo+1);CHKERRQ(ierr); 1456 ierr = PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 1457 } 1458 if (!tjsch->rctx2->reverseonestep && tjsch->rctx2->stepsleft > 0) tjsch->rctx2->stepsleft--; 1459 } 1460 /* fill stack */ 1461 if (stack->solution_only) { 1462 if (tjsch->save_stack) { 1463 if (restoredstridenum < stridenum) { 1464 ierr = StackLoadLast(tj,ts,stack,restoredstridenum);CHKERRQ(ierr); 1465 } else { 1466 ierr = StackLoadAll(tj,ts,stack,restoredstridenum);CHKERRQ(ierr); 1467 } 1468 /* recompute one step ahead */ 1469 tjsch->skip_trajectory = PETSC_TRUE; 1470 ierr = TurnForward(ts);CHKERRQ(ierr); 1471 ierr = ReCompute(ts,tjsch,stridenum*tjsch->stride-1,stridenum*tjsch->stride);CHKERRQ(ierr); 1472 tjsch->skip_trajectory = PETSC_FALSE; 1473 if (restoredstridenum < stridenum) { 1474 ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr); 1475 ierr = TurnForward(ts);CHKERRQ(ierr); 1476 ierr = ReCompute(ts,tjsch,restoredstridenum*tjsch->stride,stepnum);CHKERRQ(ierr); 1477 } else { /* stack ready, fast forward revolve status */ 1478 ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr); 1479 ierr = FastForwardRevolve(tjsch->rctx);CHKERRQ(ierr); 1480 } 1481 } else { 1482 ierr = LoadSingle(tj,ts,stack,restoredstridenum);CHKERRQ(ierr); 1483 ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr); 1484 ierr = TurnForward(ts);CHKERRQ(ierr); 1485 ierr = ReCompute(ts,tjsch,(restoredstridenum-1)*tjsch->stride,stepnum);CHKERRQ(ierr); 1486 } 1487 } else { 1488 if (tjsch->save_stack) { 1489 if (restoredstridenum < stridenum) { 1490 ierr = StackLoadLast(tj,ts,stack,restoredstridenum);CHKERRQ(ierr); 1491 /* reset revolve */ 1492 ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr); 1493 ierr = TurnForward(ts);CHKERRQ(ierr); 1494 ierr = ReCompute(ts,tjsch,restoredstridenum*tjsch->stride,stepnum);CHKERRQ(ierr); 1495 } else { /* stack ready, fast forward revolve status */ 1496 ierr = StackLoadAll(tj,ts,stack,restoredstridenum);CHKERRQ(ierr); 1497 ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr); 1498 ierr = FastForwardRevolve(tjsch->rctx);CHKERRQ(ierr); 1499 } 1500 } else { 1501 ierr = LoadSingle(tj,ts,stack,restoredstridenum);CHKERRQ(ierr); 1502 ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr); 1503 /* push first element to stack */ 1504 if (tjsch->store_stride || tjsch->rctx2->reverseonestep) { 1505 CheckpointType cptype = SOLUTION_STAGES; 1506 shift = (restoredstridenum-1)*tjsch->stride-localstepnum; 1507 ierr = PetscRevolveIntCast(tjsch->total_steps,&rtotal_steps);CHKERRQ(ierr); 1508 ierr = PetscRevolveIntCast((restoredstridenum-1)*tjsch->stride+1,&rstepnum);CHKERRQ(ierr); 1509 ierr = ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,rtotal_steps,rstepnum,1,PETSC_FALSE,&store);CHKERRQ(ierr); 1510 if (tj->monitor) { 1511 ierr = PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 1512 ierr = PetscViewerASCIIPrintf(tj->monitor,"Skip the step from %D to %D (stage values already checkpointed)\n",(restoredstridenum-1)*tjsch->stride,(restoredstridenum-1)*tjsch->stride+1);CHKERRQ(ierr); 1513 ierr = PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 1514 } 1515 ierr = ElementCreate(ts,cptype,stack,&e);CHKERRQ(ierr); 1516 ierr = ElementSet(ts,stack,&e,(restoredstridenum-1)*tjsch->stride+1,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 1517 ierr = StackPush(stack,e);CHKERRQ(ierr); 1518 } 1519 ierr = TurnForward(ts);CHKERRQ(ierr); 1520 ierr = ReCompute(ts,tjsch,(restoredstridenum-1)*tjsch->stride+1,stepnum);CHKERRQ(ierr); 1521 } 1522 } 1523 if (restoredstridenum == stridenum) diskstack->top--; 1524 tjsch->rctx->reverseonestep = PETSC_FALSE; 1525 PetscFunctionReturn(0); 1526 } 1527 1528 if (stack->solution_only) { 1529 /* restore a checkpoint */ 1530 ierr = StackTop(stack,&e);CHKERRQ(ierr); 1531 ierr = UpdateTS(ts,stack,e,stepnum,PETSC_TRUE);CHKERRQ(ierr); 1532 /* start with restoring a checkpoint */ 1533 ierr = PetscRevolveIntCast(stepnum,&rstepnum);CHKERRQ(ierr); 1534 ierr = PetscRevolveIntCast(localstepnum,&rlocalstepnum);CHKERRQ(ierr); 1535 tjsch->rctx->capo = rstepnum; 1536 tjsch->rctx->oldcapo = tjsch->rctx->capo; 1537 shift = rstepnum-rlocalstepnum; 1538 whattodo = revolve_action(&tjsch->rctx->check,&tjsch->rctx->capo,&tjsch->rctx->fine,tjsch->rctx->snaps_in,&tjsch->rctx->info,&tjsch->rctx->where); 1539 ierr = printwhattodo(tj->monitor,whattodo,tjsch->rctx,shift);CHKERRQ(ierr); 1540 ierr = TurnForward(ts);CHKERRQ(ierr); 1541 ierr = ReCompute(ts,tjsch,e->stepnum,stepnum);CHKERRQ(ierr); 1542 if (e->stepnum+1 == stepnum) { 1543 ierr = StackPop(stack,&e);CHKERRQ(ierr); 1544 } 1545 } else { 1546 PetscRevolveInt rlocalstepnum; 1547 /* restore a checkpoint */ 1548 ierr = StackTop(stack,&e);CHKERRQ(ierr); 1549 ierr = UpdateTS(ts,stack,e,stepnum,PETSC_TRUE);CHKERRQ(ierr); 1550 /* 2 revolve actions: restore a checkpoint and then advance */ 1551 ierr = PetscRevolveIntCast(tjsch->total_steps,&rtotal_steps);CHKERRQ(ierr); 1552 ierr = PetscRevolveIntCast(stridenum,&rstepnum);CHKERRQ(ierr); 1553 ierr = PetscRevolveIntCast(localstepnum,&rlocalstepnum);CHKERRQ(ierr); 1554 ierr = ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,rtotal_steps,rstepnum,rlocalstepnum,PETSC_FALSE,&store);CHKERRQ(ierr); 1555 if (tj->monitor) { 1556 ierr = PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 1557 ierr = 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);CHKERRQ(ierr); 1558 ierr = PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 1559 } 1560 if (!tjsch->rctx->reverseonestep && tjsch->rctx->stepsleft > 0) tjsch->rctx->stepsleft--; 1561 if (e->stepnum < stepnum) { 1562 ierr = TurnForward(ts);CHKERRQ(ierr); 1563 ierr = ReCompute(ts,tjsch,e->stepnum,stepnum);CHKERRQ(ierr); 1564 } 1565 if (e->stepnum == stepnum) { 1566 ierr = StackPop(stack,&e);CHKERRQ(ierr); 1567 } 1568 } 1569 tjsch->rctx->reverseonestep = PETSC_FALSE; 1570 PetscFunctionReturn(0); 1571 } 1572 1573 static PetscErrorCode SetTrajRMS(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X) 1574 { 1575 Stack *stack = &tjsch->stack; 1576 PetscInt store; 1577 StackElement e; 1578 PetscRevolveInt rtotal_steps,rstepnum; 1579 PetscErrorCode ierr; 1580 1581 PetscFunctionBegin; 1582 if (!stack->solution_only && stepnum == 0) PetscFunctionReturn(0); 1583 if (stack->solution_only && stepnum == tjsch->total_steps) PetscFunctionReturn(0); 1584 ierr = PetscRevolveIntCast(tjsch->total_steps,&rtotal_steps);CHKERRQ(ierr); 1585 ierr = PetscRevolveIntCast(stepnum,&rstepnum);CHKERRQ(ierr); 1586 ierr = ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,rtotal_steps,rstepnum,rstepnum,PETSC_FALSE,&store);CHKERRQ(ierr); 1587 if (store == 1) { 1588 CheckpointType cptype; 1589 if (stepnum < stack->top) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_MEMC,"Illegal modification of a non-top stack element"); 1590 cptype = stack->solution_only ? SOLUTIONONLY : SOLUTION_STAGES; 1591 ierr = ElementCreate(ts,cptype,stack,&e);CHKERRQ(ierr); 1592 ierr = ElementSet(ts,stack,&e,stepnum,time,X);CHKERRQ(ierr); 1593 ierr = StackPush(stack,e);CHKERRQ(ierr); 1594 } else if (store == 2) { 1595 ierr = DumpSingle(tj,ts,stack,tjsch->rctx->check+1);CHKERRQ(ierr); 1596 } 1597 PetscFunctionReturn(0); 1598 } 1599 1600 static PetscErrorCode GetTrajRMS(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum) 1601 { 1602 Stack *stack = &tjsch->stack; 1603 PetscRevolveInt whattodo,shift,rstepnum; 1604 PetscInt restart; 1605 PetscBool ondisk; 1606 StackElement e; 1607 PetscErrorCode ierr; 1608 1609 PetscFunctionBegin; 1610 if (stepnum == 0 || stepnum == tjsch->total_steps) { 1611 ierr = TurnBackward(ts);CHKERRQ(ierr); 1612 tjsch->rctx->reverseonestep = PETSC_FALSE; 1613 PetscFunctionReturn(0); 1614 } 1615 ierr = PetscRevolveIntCast(stepnum,&rstepnum);CHKERRQ(ierr); 1616 tjsch->rctx->capo = rstepnum; 1617 tjsch->rctx->oldcapo = tjsch->rctx->capo; 1618 shift = 0; 1619 whattodo = revolve_action(&tjsch->rctx->check,&tjsch->rctx->capo,&tjsch->rctx->fine,tjsch->rctx->snaps_in,&tjsch->rctx->info,&tjsch->rctx->where); /* whattodo=restore */ 1620 ierr = printwhattodo(tj->monitor,whattodo,tjsch->rctx,shift);CHKERRQ(ierr); 1621 /* restore a checkpoint */ 1622 restart = tjsch->rctx->capo; 1623 if (!tjsch->rctx->where) { 1624 ondisk = PETSC_TRUE; 1625 ierr = LoadSingle(tj,ts,stack,tjsch->rctx->check+1);CHKERRQ(ierr); 1626 ierr = TurnBackward(ts);CHKERRQ(ierr); 1627 } else { 1628 ondisk = PETSC_FALSE; 1629 ierr = StackTop(stack,&e);CHKERRQ(ierr); 1630 ierr = UpdateTS(ts,stack,e,stepnum,PETSC_TRUE);CHKERRQ(ierr); 1631 } 1632 if (!stack->solution_only) { /* whattodo must be 5 or 8 */ 1633 /* ask Revolve what to do next */ 1634 tjsch->rctx->oldcapo = tjsch->rctx->capo; 1635 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*/ 1636 ierr = printwhattodo(tj->monitor,whattodo,tjsch->rctx,shift);CHKERRQ(ierr); 1637 if (whattodo == 3 || whattodo == 4) tjsch->rctx->reverseonestep = PETSC_TRUE; 1638 if (whattodo == 1) tjsch->rctx->stepsleft = tjsch->rctx->capo-tjsch->rctx->oldcapo; 1639 if (tj->monitor) { 1640 ierr = PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 1641 ierr = PetscViewerASCIIPrintf(tj->monitor,"Skip the step from %D to %D (stage values already checkpointed)\n",tjsch->rctx->oldcapo,tjsch->rctx->oldcapo+1);CHKERRQ(ierr); 1642 ierr = PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 1643 } 1644 if (!tjsch->rctx->reverseonestep && tjsch->rctx->stepsleft > 0) tjsch->rctx->stepsleft--; 1645 restart++; /* skip one step */ 1646 } 1647 if (stack->solution_only || (!stack->solution_only && restart < stepnum)) { 1648 ierr = TurnForward(ts);CHKERRQ(ierr); 1649 ierr = ReCompute(ts,tjsch,restart,stepnum);CHKERRQ(ierr); 1650 } 1651 if (!ondisk && ( (stack->solution_only && e->stepnum+1 == stepnum) || (!stack->solution_only && e->stepnum == stepnum))) { 1652 ierr = StackPop(stack,&e);CHKERRQ(ierr); 1653 } 1654 tjsch->rctx->reverseonestep = PETSC_FALSE; 1655 PetscFunctionReturn(0); 1656 } 1657 #endif 1658 1659 #if defined(PETSC_HAVE_CAMS) 1660 /* Optimal offline adjoint checkpointing for multistage time integration methods */ 1661 static PetscErrorCode SetTrajAOF(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X) 1662 { 1663 Stack *stack = &tjsch->stack; 1664 StackElement e; 1665 PetscErrorCode ierr; 1666 1667 PetscFunctionBegin; 1668 /* skip if no checkpoint to use. This also avoids an error when num_units_avail=0 */ 1669 if (tjsch->actx->nextcheckpointstep == -1) PetscFunctionReturn(0); 1670 if (stepnum == 0) { /* When placing the first checkpoint, no need to change the units available */ 1671 if (stack->solution_only) { 1672 ierr = offline_ca(tjsch->actx->lastcheckpointstep,tjsch->actx->num_units_avail,tjsch->actx->endstep,&tjsch->actx->nextcheckpointstep);CHKERRQ(ierr); 1673 } else { 1674 /* First two arguments must be -1 when first time calling cams */ 1675 ierr = offline_cams(tjsch->actx->lastcheckpointstep,tjsch->actx->lastcheckpointtype,tjsch->actx->num_units_avail,tjsch->actx->endstep,tjsch->actx->num_stages,&tjsch->actx->nextcheckpointstep,&tjsch->actx->nextcheckpointtype);CHKERRQ(ierr); 1676 } 1677 } 1678 1679 if (stack->solution_only && stepnum == tjsch->total_steps) PetscFunctionReturn(0); 1680 1681 if (tjsch->actx->nextcheckpointstep == stepnum) { 1682 if (stepnum < stack->top) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_MEMC,"Illegal modification of a non-top stack element"); 1683 1684 if (tjsch->actx->nextcheckpointtype == 2) { /* solution + stage values */ 1685 if (tj->monitor) { 1686 ierr = PetscViewerASCIIPrintf(tj->monitor,"Store in checkpoint number %D with stage values and solution (located in RAM)\n",stepnum);CHKERRQ(ierr); 1687 } 1688 ierr = ElementCreate(ts,SOLUTION_STAGES,stack,&e);CHKERRQ(ierr); 1689 ierr = ElementSet(ts,stack,&e,stepnum,time,X);CHKERRQ(ierr); 1690 } 1691 if (tjsch->actx->nextcheckpointtype == 1) { 1692 if (tj->monitor) { 1693 ierr = PetscViewerASCIIPrintf(tj->monitor,"Store in checkpoint number %D with stage values (located in RAM)\n",stepnum);CHKERRQ(ierr); 1694 } 1695 ierr = ElementCreate(ts,STAGESONLY,stack,&e);CHKERRQ(ierr); 1696 ierr = ElementSet(ts,stack,&e,stepnum,time,X);CHKERRQ(ierr); 1697 } 1698 if (tjsch->actx->nextcheckpointtype == 0) { /* solution only */ 1699 if (tj->monitor) { 1700 ierr = PetscViewerASCIIPrintf(tj->monitor,"Store in checkpoint number %D (located in RAM)\n",stepnum);CHKERRQ(ierr); 1701 } 1702 ierr = ElementCreate(ts,SOLUTIONONLY,stack,&e);CHKERRQ(ierr); 1703 ierr = ElementSet(ts,stack,&e,stepnum,time,X);CHKERRQ(ierr); 1704 } 1705 ierr = StackPush(stack,e);CHKERRQ(ierr); 1706 1707 tjsch->actx->lastcheckpointstep = stepnum; 1708 if (stack->solution_only) { 1709 ierr = offline_ca(tjsch->actx->lastcheckpointstep,tjsch->actx->num_units_avail,tjsch->actx->endstep,&tjsch->actx->nextcheckpointstep);CHKERRQ(ierr); 1710 tjsch->actx->num_units_avail--; 1711 } else { 1712 tjsch->actx->lastcheckpointtype = tjsch->actx->nextcheckpointtype; 1713 ierr = offline_cams(tjsch->actx->lastcheckpointstep,tjsch->actx->lastcheckpointtype,tjsch->actx->num_units_avail,tjsch->actx->endstep,tjsch->actx->num_stages,&tjsch->actx->nextcheckpointstep,&tjsch->actx->nextcheckpointtype);CHKERRQ(ierr); 1714 if (tjsch->actx->lastcheckpointtype == 2) tjsch->actx->num_units_avail -= tjsch->actx->num_stages+1; 1715 if (tjsch->actx->lastcheckpointtype == 1) tjsch->actx->num_units_avail -= tjsch->actx->num_stages; 1716 if (tjsch->actx->lastcheckpointtype == 0) tjsch->actx->num_units_avail--; 1717 } 1718 } 1719 PetscFunctionReturn(0); 1720 } 1721 1722 static PetscErrorCode GetTrajAOF(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum) 1723 { 1724 Stack *stack = &tjsch->stack; 1725 StackElement e; 1726 PetscInt estepnum; 1727 PetscErrorCode ierr; 1728 1729 PetscFunctionBegin; 1730 if (stepnum == 0 || stepnum == tjsch->total_steps) { 1731 ierr = TurnBackward(ts);CHKERRQ(ierr); 1732 PetscFunctionReturn(0); 1733 } 1734 /* Restore a checkpoint */ 1735 ierr = StackTop(stack,&e);CHKERRQ(ierr); 1736 estepnum = e->stepnum; 1737 if (estepnum == stepnum && e->cptype == SOLUTIONONLY) { /* discard the checkpoint if not useful (corner case) */ 1738 ierr = StackPop(stack,&e);CHKERRQ(ierr); 1739 tjsch->actx->num_units_avail++; 1740 ierr = StackTop(stack,&e);CHKERRQ(ierr); 1741 estepnum = e->stepnum; 1742 } 1743 /* Update TS with stage values if an adjoint step can be taken immediately */ 1744 if (HaveStages(e->cptype)) { 1745 if (tj->monitor) { 1746 ierr = PetscViewerASCIIPrintf(tj->monitor,"Restore in checkpoint number %D with stage values (located in RAM)\n",e->stepnum);CHKERRQ(ierr); 1747 } 1748 if (e->cptype == STAGESONLY) tjsch->actx->num_units_avail += tjsch->actx->num_stages; 1749 if (e->cptype == SOLUTION_STAGES) tjsch->actx->num_units_avail += tjsch->actx->num_stages+1; 1750 } else { 1751 if (tj->monitor) { 1752 ierr = PetscViewerASCIIPrintf(tj->monitor,"Restore in checkpoint number %D (located in RAM)\n",e->stepnum);CHKERRQ(ierr); 1753 } 1754 tjsch->actx->num_units_avail++; 1755 } 1756 ierr = UpdateTS(ts,stack,e,stepnum,PETSC_TRUE);CHKERRQ(ierr); 1757 /* Query the scheduler */ 1758 tjsch->actx->lastcheckpointstep = estepnum; 1759 tjsch->actx->endstep = stepnum; 1760 if (stack->solution_only) { /* start with restoring a checkpoint */ 1761 ierr = offline_ca(tjsch->actx->lastcheckpointstep,tjsch->actx->num_units_avail,tjsch->actx->endstep,&tjsch->actx->nextcheckpointstep);CHKERRQ(ierr); 1762 } else { /* 2 revolve actions: restore a checkpoint and then advance */ 1763 tjsch->actx->lastcheckpointtype = e->cptype; 1764 ierr = offline_cams(tjsch->actx->lastcheckpointstep,tjsch->actx->lastcheckpointtype,tjsch->actx->num_units_avail,tjsch->actx->endstep,tjsch->actx->num_stages,&tjsch->actx->nextcheckpointstep, &tjsch->actx->nextcheckpointtype);CHKERRQ(ierr); 1765 } 1766 /* Discard the checkpoint if not needed, decrease the number of available checkpoints if it still stays in stack */ 1767 if (HaveStages(e->cptype)) { 1768 if (estepnum == stepnum) { 1769 ierr = StackPop(stack,&e);CHKERRQ(ierr); 1770 } else { 1771 if (e->cptype == STAGESONLY) tjsch->actx->num_units_avail -= tjsch->actx->num_stages; 1772 if (e->cptype == SOLUTION_STAGES) tjsch->actx->num_units_avail -= tjsch->actx->num_stages+1; 1773 } 1774 } else { 1775 if (estepnum+1 == stepnum) { 1776 ierr = StackPop(stack,&e);CHKERRQ(ierr); 1777 } else { 1778 tjsch->actx->num_units_avail--; 1779 } 1780 } 1781 /* Recompute from the restored checkpoint */ 1782 if (stack->solution_only || (!stack->solution_only && estepnum < stepnum)) { 1783 ierr = TurnForward(ts);CHKERRQ(ierr); 1784 ierr = ReCompute(ts,tjsch,estepnum,stepnum);CHKERRQ(ierr); 1785 } 1786 PetscFunctionReturn(0); 1787 } 1788 #endif 1789 1790 static PetscErrorCode TSTrajectorySet_Memory(TSTrajectory tj,TS ts,PetscInt stepnum,PetscReal time,Vec X) 1791 { 1792 TJScheduler *tjsch = (TJScheduler*)tj->data; 1793 PetscErrorCode ierr; 1794 1795 PetscFunctionBegin; 1796 if (!tjsch->recompute) { /* use global stepnum in the forward sweep */ 1797 ierr = TSGetStepNumber(ts,&stepnum);CHKERRQ(ierr); 1798 } 1799 /* for consistency */ 1800 if (!tjsch->recompute && stepnum == 0) ts->ptime_prev = ts->ptime-ts->time_step; 1801 switch (tjsch->stype) { 1802 case NONE: 1803 if (tj->adjoint_solve_mode) { 1804 ierr = SetTrajN(ts,tjsch,stepnum,time,X);CHKERRQ(ierr); 1805 } else { 1806 ierr = SetTrajN_2(ts,tjsch,stepnum,time,X);CHKERRQ(ierr); 1807 } 1808 break; 1809 case TWO_LEVEL_NOREVOLVE: 1810 if (!tj->adjoint_solve_mode) SETERRQ(PetscObjectComm((PetscObject)tj),PETSC_ERR_SUP,"Not implemented"); 1811 ierr = SetTrajTLNR(tj,ts,tjsch,stepnum,time,X);CHKERRQ(ierr); 1812 break; 1813 #if defined(PETSC_HAVE_REVOLVE) 1814 case TWO_LEVEL_REVOLVE: 1815 if (!tj->adjoint_solve_mode) SETERRQ(PetscObjectComm((PetscObject)tj),PETSC_ERR_SUP,"Not implemented"); 1816 ierr = SetTrajTLR(tj,ts,tjsch,stepnum,time,X);CHKERRQ(ierr); 1817 break; 1818 case TWO_LEVEL_TWO_REVOLVE: 1819 if (!tj->adjoint_solve_mode) SETERRQ(PetscObjectComm((PetscObject)tj),PETSC_ERR_SUP,"Not implemented"); 1820 ierr = SetTrajTLTR(tj,ts,tjsch,stepnum,time,X);CHKERRQ(ierr); 1821 break; 1822 case REVOLVE_OFFLINE: 1823 if (!tj->adjoint_solve_mode) SETERRQ(PetscObjectComm((PetscObject)tj),PETSC_ERR_SUP,"Not implemented"); 1824 ierr = SetTrajROF(tj,ts,tjsch,stepnum,time,X);CHKERRQ(ierr); 1825 break; 1826 case REVOLVE_ONLINE: 1827 if (!tj->adjoint_solve_mode) SETERRQ(PetscObjectComm((PetscObject)tj),PETSC_ERR_SUP,"Not implemented"); 1828 ierr = SetTrajRON(tj,ts,tjsch,stepnum,time,X);CHKERRQ(ierr); 1829 break; 1830 case REVOLVE_MULTISTAGE: 1831 if (!tj->adjoint_solve_mode) SETERRQ(PetscObjectComm((PetscObject)tj),PETSC_ERR_SUP,"Not implemented"); 1832 ierr = SetTrajRMS(tj,ts,tjsch,stepnum,time,X);CHKERRQ(ierr); 1833 break; 1834 #endif 1835 #if defined(PETSC_HAVE_CAMS) 1836 case CAMS_OFFLINE: 1837 if (!tj->adjoint_solve_mode) SETERRQ(PetscObjectComm((PetscObject)tj),PETSC_ERR_SUP,"Not implemented"); 1838 ierr = SetTrajAOF(tj,ts,tjsch,stepnum,time,X);CHKERRQ(ierr); 1839 break; 1840 #endif 1841 default: 1842 break; 1843 } 1844 PetscFunctionReturn(0); 1845 } 1846 1847 static PetscErrorCode TSTrajectoryGet_Memory(TSTrajectory tj,TS ts,PetscInt stepnum,PetscReal *t) 1848 { 1849 TJScheduler *tjsch = (TJScheduler*)tj->data; 1850 PetscErrorCode ierr; 1851 1852 PetscFunctionBegin; 1853 if (tj->adjoint_solve_mode && stepnum == 0) { 1854 ierr = TSTrajectoryReset(tj);CHKERRQ(ierr); /* reset TSTrajectory so users do not need to reset TSTrajectory */ 1855 PetscFunctionReturn(0); 1856 } 1857 switch (tjsch->stype) { 1858 case NONE: 1859 if (tj->adjoint_solve_mode) { 1860 ierr = GetTrajN(ts,tjsch,stepnum);CHKERRQ(ierr); 1861 } else { 1862 ierr = GetTrajN_2(ts,tjsch,stepnum);CHKERRQ(ierr); 1863 } 1864 break; 1865 case TWO_LEVEL_NOREVOLVE: 1866 if (!tj->adjoint_solve_mode) SETERRQ(PetscObjectComm((PetscObject)tj),PETSC_ERR_SUP,"Not implemented"); 1867 ierr = GetTrajTLNR(tj,ts,tjsch,stepnum);CHKERRQ(ierr); 1868 break; 1869 #if defined(PETSC_HAVE_REVOLVE) 1870 case TWO_LEVEL_REVOLVE: 1871 if (!tj->adjoint_solve_mode) SETERRQ(PetscObjectComm((PetscObject)tj),PETSC_ERR_SUP,"Not implemented"); 1872 ierr = GetTrajTLR(tj,ts,tjsch,stepnum);CHKERRQ(ierr); 1873 break; 1874 case TWO_LEVEL_TWO_REVOLVE: 1875 if (!tj->adjoint_solve_mode) SETERRQ(PetscObjectComm((PetscObject)tj),PETSC_ERR_SUP,"Not implemented"); 1876 ierr = GetTrajTLTR(tj,ts,tjsch,stepnum);CHKERRQ(ierr); 1877 break; 1878 case REVOLVE_OFFLINE: 1879 if (!tj->adjoint_solve_mode) SETERRQ(PetscObjectComm((PetscObject)tj),PETSC_ERR_SUP,"Not implemented"); 1880 ierr = GetTrajROF(tj,ts,tjsch,stepnum);CHKERRQ(ierr); 1881 break; 1882 case REVOLVE_ONLINE: 1883 if (!tj->adjoint_solve_mode) SETERRQ(PetscObjectComm((PetscObject)tj),PETSC_ERR_SUP,"Not implemented"); 1884 ierr = GetTrajRON(tj,ts,tjsch,stepnum);CHKERRQ(ierr); 1885 break; 1886 case REVOLVE_MULTISTAGE: 1887 if (!tj->adjoint_solve_mode) SETERRQ(PetscObjectComm((PetscObject)tj),PETSC_ERR_SUP,"Not implemented"); 1888 ierr = GetTrajRMS(tj,ts,tjsch,stepnum);CHKERRQ(ierr); 1889 break; 1890 #endif 1891 #if defined(PETSC_HAVE_CAMS) 1892 case CAMS_OFFLINE: 1893 if (!tj->adjoint_solve_mode) SETERRQ(PetscObjectComm((PetscObject)tj),PETSC_ERR_SUP,"Not implemented"); 1894 ierr = GetTrajAOF(tj,ts,tjsch,stepnum);CHKERRQ(ierr); 1895 break; 1896 #endif 1897 default: 1898 break; 1899 } 1900 PetscFunctionReturn(0); 1901 } 1902 1903 PETSC_UNUSED static PetscErrorCode TSTrajectorySetStride_Memory(TSTrajectory tj,PetscInt stride) 1904 { 1905 TJScheduler *tjsch = (TJScheduler*)tj->data; 1906 1907 PetscFunctionBegin; 1908 tjsch->stride = stride; 1909 PetscFunctionReturn(0); 1910 } 1911 1912 static PetscErrorCode TSTrajectorySetMaxCpsRAM_Memory(TSTrajectory tj,PetscInt max_cps_ram) 1913 { 1914 TJScheduler *tjsch = (TJScheduler*)tj->data; 1915 1916 PetscFunctionBegin; 1917 tjsch->max_cps_ram = max_cps_ram; 1918 PetscFunctionReturn(0); 1919 } 1920 1921 static PetscErrorCode TSTrajectorySetMaxCpsDisk_Memory(TSTrajectory tj,PetscInt max_cps_disk) 1922 { 1923 TJScheduler *tjsch = (TJScheduler*)tj->data; 1924 1925 PetscFunctionBegin; 1926 tjsch->max_cps_disk = max_cps_disk; 1927 PetscFunctionReturn(0); 1928 } 1929 1930 static PetscErrorCode TSTrajectorySetMaxUnitsRAM_Memory(TSTrajectory tj,PetscInt max_units_ram) 1931 { 1932 TJScheduler *tjsch = (TJScheduler*)tj->data; 1933 1934 PetscFunctionBegin; 1935 if (!tjsch->max_cps_ram) SETERRQ(PetscObjectComm((PetscObject)tj),PETSC_ERR_ARG_INCOMP,"Conflict with -ts_trjaectory_max_cps_ram or TSTrajectorySetMaxCpsRAM. You can set max_cps_ram or max_units_ram, but not both at the same time."); 1936 tjsch->max_units_ram = max_units_ram; 1937 PetscFunctionReturn(0); 1938 } 1939 1940 static PetscErrorCode TSTrajectorySetMaxUnitsDisk_Memory(TSTrajectory tj,PetscInt max_units_disk) 1941 { 1942 TJScheduler *tjsch = (TJScheduler*)tj->data; 1943 1944 PetscFunctionBegin; 1945 if (!tjsch->max_cps_disk) SETERRQ(PetscObjectComm((PetscObject)tj),PETSC_ERR_ARG_INCOMP,"Conflict with -ts_trjaectory_max_cps_disk or TSTrajectorySetMaxCpsDisk. You can set max_cps_disk or max_units_disk, but not both at the same time."); 1946 tjsch->max_units_ram = max_units_disk; 1947 PetscFunctionReturn(0); 1948 } 1949 1950 #if defined(PETSC_HAVE_REVOLVE) 1951 PETSC_UNUSED static PetscErrorCode TSTrajectorySetRevolveOnline(TSTrajectory tj,PetscBool use_online) 1952 { 1953 TJScheduler *tjsch = (TJScheduler*)tj->data; 1954 1955 PetscFunctionBegin; 1956 tjsch->use_online = use_online; 1957 PetscFunctionReturn(0); 1958 } 1959 #endif 1960 1961 PETSC_UNUSED static PetscErrorCode TSTrajectorySetSaveStack(TSTrajectory tj,PetscBool save_stack) 1962 { 1963 TJScheduler *tjsch = (TJScheduler*)tj->data; 1964 1965 PetscFunctionBegin; 1966 tjsch->save_stack = save_stack; 1967 PetscFunctionReturn(0); 1968 } 1969 1970 PETSC_UNUSED static PetscErrorCode TSTrajectorySetUseDRAM(TSTrajectory tj,PetscBool use_dram) 1971 { 1972 TJScheduler *tjsch = (TJScheduler*)tj->data; 1973 1974 PetscFunctionBegin; 1975 tjsch->stack.use_dram = use_dram; 1976 PetscFunctionReturn(0); 1977 } 1978 1979 /*@C 1980 TSTrajectorySetMaxCpsRAM - Set maximum number of checkpoints in RAM 1981 1982 Logically collective 1983 1984 Input Parameter: 1985 . tj - tstrajectory context 1986 1987 Output Parameter: 1988 . max_cps_ram - maximum number of checkpoints in RAM 1989 1990 Level: intermediate 1991 1992 .seealso: TSTrajectorySetMaxUnitsRAM() 1993 @*/ 1994 PetscErrorCode TSTrajectorySetMaxCpsRAM(TSTrajectory tj,PetscInt max_cps_ram) 1995 { 1996 PetscErrorCode ierr; 1997 1998 PetscFunctionBegin; 1999 ierr = PetscUseMethod(tj,"TSTrajectorySetMaxCpsRAM_C",(TSTrajectory,PetscInt),(tj,max_cps_ram));CHKERRQ(ierr); 2000 PetscFunctionReturn(0); 2001 } 2002 2003 /*@C 2004 TSTrajectorySetMaxCpsDisk - Set maximum number of checkpoints on disk 2005 2006 Logically collective 2007 2008 Input Parameter: 2009 . tj - tstrajectory context 2010 2011 Output Parameter: 2012 . max_cps_disk - maximum number of checkpoints on disk 2013 2014 Level: intermediate 2015 2016 .seealso: TSTrajectorySetMaxUnitsDisk(), TSTrajectorySetMaxUnitsRAM() 2017 @*/ 2018 PetscErrorCode TSTrajectorySetMaxCpsDisk(TSTrajectory tj,PetscInt max_cps_disk) 2019 { 2020 PetscErrorCode ierr; 2021 2022 PetscFunctionBegin; 2023 ierr = PetscUseMethod(tj,"TSTrajectorySetMaxCpsDisk_C",(TSTrajectory,PetscInt),(tj,max_cps_disk));CHKERRQ(ierr); 2024 PetscFunctionReturn(0); 2025 } 2026 2027 /*@C 2028 TSTrajectorySetMaxUnitsRAM - Set maximum number of checkpointing units in RAM 2029 2030 Logically collective 2031 2032 Input Parameter: 2033 . tj - tstrajectory context 2034 2035 Output Parameter: 2036 . max_units_ram - maximum number of checkpointing units in RAM 2037 2038 Level: intermediate 2039 2040 .seealso: TSTrajectorySetMaxCpsRAM() 2041 @*/ 2042 PetscErrorCode TSTrajectorySetMaxUnitsRAM(TSTrajectory tj,PetscInt max_units_ram) 2043 { 2044 PetscErrorCode ierr; 2045 2046 PetscFunctionBegin; 2047 ierr = PetscUseMethod(tj,"TSTrajectorySetMaxUnitsRAM_C",(TSTrajectory,PetscInt),(tj,max_units_ram));CHKERRQ(ierr); 2048 PetscFunctionReturn(0); 2049 } 2050 2051 /*@C 2052 TSTrajectorySetMaxUnitsDisk - Set maximum number of checkpointing units on disk 2053 2054 Logically collective 2055 2056 Input Parameter: 2057 . tj - tstrajectory context 2058 2059 Output Parameter: 2060 . max_units_disk - maximum number of checkpointing units on disk 2061 2062 Level: intermediate 2063 2064 .seealso: TSTrajectorySetMaxCpsDisk() 2065 @*/ 2066 PetscErrorCode TSTrajectorySetMaxUnitsDisk(TSTrajectory tj,PetscInt max_units_disk) 2067 { 2068 PetscErrorCode ierr; 2069 2070 PetscFunctionBegin; 2071 ierr = PetscUseMethod(tj,"TSTrajectorySetMaxUnitsDisk_C",(TSTrajectory,PetscInt),(tj,max_units_disk));CHKERRQ(ierr); 2072 PetscFunctionReturn(0); 2073 } 2074 2075 static PetscErrorCode TSTrajectorySetFromOptions_Memory(PetscOptionItems *PetscOptionsObject,TSTrajectory tj) 2076 { 2077 TJScheduler *tjsch = (TJScheduler*)tj->data; 2078 PetscInt max_cps_ram,max_cps_disk,max_units_ram,max_units_disk; 2079 PetscBool flg; 2080 PetscErrorCode ierr; 2081 2082 PetscFunctionBegin; 2083 ierr = PetscOptionsHead(PetscOptionsObject,"Memory based TS trajectory options");CHKERRQ(ierr); 2084 { 2085 ierr = PetscOptionsInt("-ts_trajectory_max_cps_ram","Maximum number of checkpoints in RAM","TSTrajectorySetMaxCpsRAM",tjsch->max_cps_ram,&max_cps_ram,&flg);CHKERRQ(ierr); 2086 if (flg) { 2087 ierr = TSTrajectorySetMaxCpsRAM(tj,max_cps_ram);CHKERRQ(ierr); 2088 } 2089 ierr = PetscOptionsInt("-ts_trajectory_max_cps_disk","Maximum number of checkpoints on disk","TSTrajectorySetMaxCpsDisk",tjsch->max_cps_disk,&max_cps_disk,&flg);CHKERRQ(ierr); 2090 if (flg) { 2091 ierr = TSTrajectorySetMaxCpsDisk(tj,max_cps_disk);CHKERRQ(ierr); 2092 } 2093 ierr = PetscOptionsInt("-ts_trajectory_max_units_ram","Maximum number of checkpointing units in RAM","TSTrajectorySetMaxUnitsRAM",tjsch->max_units_ram,&max_units_ram,&flg);CHKERRQ(ierr); 2094 if (flg) { 2095 ierr = TSTrajectorySetMaxUnitsRAM(tj,max_units_ram);CHKERRQ(ierr); 2096 } 2097 ierr = PetscOptionsInt("-ts_trajectory_max_units_disk","Maximum number of checkpointing units on disk","TSTrajectorySetMaxUnitsDisk",tjsch->max_units_disk,&max_units_disk,&flg);CHKERRQ(ierr); 2098 if (flg) { 2099 ierr = TSTrajectorySetMaxUnitsDisk(tj,max_units_disk);CHKERRQ(ierr); 2100 } 2101 ierr = PetscOptionsInt("-ts_trajectory_stride","Stride to save checkpoints to file","TSTrajectorySetStride",tjsch->stride,&tjsch->stride,NULL);CHKERRQ(ierr); 2102 #if defined(PETSC_HAVE_REVOLVE) 2103 ierr = PetscOptionsBool("-ts_trajectory_revolve_online","Trick TS trajectory into using online mode of revolve","TSTrajectorySetRevolveOnline",tjsch->use_online,&tjsch->use_online,NULL);CHKERRQ(ierr); 2104 #endif 2105 ierr = PetscOptionsBool("-ts_trajectory_save_stack","Save all stack to disk","TSTrajectorySetSaveStack",tjsch->save_stack,&tjsch->save_stack,NULL);CHKERRQ(ierr); 2106 ierr = PetscOptionsBool("-ts_trajectory_use_dram","Use DRAM for checkpointing","TSTrajectorySetUseDRAM",tjsch->stack.use_dram,&tjsch->stack.use_dram,NULL);CHKERRQ(ierr); 2107 } 2108 ierr = PetscOptionsTail();CHKERRQ(ierr); 2109 tjsch->stack.solution_only = tj->solution_only; 2110 PetscFunctionReturn(0); 2111 } 2112 2113 static PetscErrorCode TSTrajectorySetUp_Memory(TSTrajectory tj,TS ts) 2114 { 2115 TJScheduler *tjsch = (TJScheduler*)tj->data; 2116 Stack *stack = &tjsch->stack; 2117 #if defined(PETSC_HAVE_REVOLVE) 2118 RevolveCTX *rctx,*rctx2; 2119 DiskStack *diskstack = &tjsch->diskstack; 2120 PetscInt diskblocks; 2121 #endif 2122 PetscInt numY,total_steps; 2123 PetscBool fixedtimestep; 2124 PetscErrorCode ierr; 2125 2126 PetscFunctionBegin; 2127 if (ts->adapt) { 2128 ierr = PetscObjectTypeCompare((PetscObject)ts->adapt,TSADAPTNONE,&fixedtimestep);CHKERRQ(ierr); 2129 } else { 2130 fixedtimestep = PETSC_TRUE; 2131 } 2132 total_steps = (PetscInt)(PetscCeilReal((ts->max_time-ts->ptime)/ts->time_step)); 2133 total_steps = total_steps < 0 ? PETSC_MAX_INT : total_steps; 2134 if (fixedtimestep) tjsch->total_steps = PetscMin(ts->max_steps,total_steps); 2135 2136 ierr = TSGetStages(ts,&numY,PETSC_IGNORE);CHKERRQ(ierr); 2137 if (stack->solution_only) { 2138 if (tjsch->max_units_ram) tjsch->max_cps_ram = tjsch->max_units_ram; 2139 else tjsch->max_units_ram = tjsch->max_cps_ram; 2140 if (tjsch->max_units_disk) tjsch->max_cps_disk = tjsch->max_units_disk; 2141 } else { 2142 if (tjsch->max_units_ram) tjsch->max_cps_ram = (ts->stifflyaccurate) ? tjsch->max_units_ram/numY : tjsch->max_units_ram/(numY+1); 2143 else tjsch->max_units_ram = (ts->stifflyaccurate) ? numY*tjsch->max_cps_ram : (numY+1)*tjsch->max_cps_ram; 2144 if (tjsch->max_units_disk) tjsch->max_cps_disk = (ts->stifflyaccurate) ? tjsch->max_units_disk/numY : tjsch->max_units_disk/(numY+1); 2145 else tjsch->max_units_disk = (ts->stifflyaccurate) ? numY*tjsch->max_cps_disk : (numY+1)*tjsch->max_cps_disk; 2146 } 2147 if (tjsch->max_cps_ram > 0) stack->stacksize = tjsch->max_units_ram; /* maximum stack size. Could be overallocated. */ 2148 2149 /* Determine the scheduler type */ 2150 if (tjsch->stride > 1) { /* two level mode */ 2151 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."); 2152 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 */ 2153 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 */ 2154 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 */ 2155 } else { /* single level mode */ 2156 if (fixedtimestep) { 2157 if (tjsch->max_cps_ram >= tjsch->total_steps-1 || tjsch->max_cps_ram == -1) 2158 tjsch->stype = NONE; /* checkpoint all */ 2159 else 2160 #if defined(PETSC_HAVE_CAMS) 2161 /* CAMS has the priority over revolve */ 2162 tjsch->stype = CAMS_OFFLINE; 2163 #else 2164 tjsch->stype = (tjsch->max_cps_disk>1) ? REVOLVE_MULTISTAGE : REVOLVE_OFFLINE; 2165 #endif 2166 } else tjsch->stype = NONE; /* checkpoint all for adaptive time step */ 2167 #if defined(PETSC_HAVE_REVOLVE) 2168 if (tjsch->use_online) tjsch->stype = REVOLVE_ONLINE; /* trick into online (for testing purpose only) */ 2169 #endif 2170 if (tjsch->stype != NONE && tjsch->max_cps_ram < 1 && tjsch->max_cps_disk < 1) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_INCOMP,"The specified storage capacity is insufficient for one checkpoint, which is the minimum"); 2171 } 2172 if (tjsch->stype >= CAMS_OFFLINE) { 2173 #ifndef PETSC_HAVE_CAMS 2174 SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"CAMS 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-cams."); 2175 #else 2176 CAMSCTX *actx; 2177 PetscInt ns = 0; 2178 if (stack->solution_only) { 2179 offline_ca_create(tjsch->total_steps,tjsch->max_cps_ram); 2180 } else { 2181 ierr = TSGetStages(ts,&ns,PETSC_IGNORE);CHKERRQ(ierr); 2182 offline_cams_create(tjsch->total_steps,tjsch->max_units_ram,ns,ts->stifflyaccurate); 2183 } 2184 ierr = PetscNew(&actx);CHKERRQ(ierr); 2185 actx->lastcheckpointstep = -1; /* -1 can trigger the initialization of CAMS */ 2186 actx->lastcheckpointtype = -1; /* -1 can trigger the initialization of CAMS */ 2187 actx->endstep = tjsch->total_steps; 2188 actx->num_units_avail = tjsch->max_units_ram; 2189 actx->num_stages = ns; 2190 tjsch->actx = actx; 2191 #endif 2192 } else if (tjsch->stype > TWO_LEVEL_NOREVOLVE) { 2193 #ifndef PETSC_HAVE_REVOLVE 2194 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."); 2195 #else 2196 PetscRevolveInt rfine,rsnaps,rsnaps2; 2197 2198 switch (tjsch->stype) { 2199 case TWO_LEVEL_REVOLVE: 2200 ierr = PetscRevolveIntCast(tjsch->stride,&rfine);CHKERRQ(ierr); 2201 ierr = PetscRevolveIntCast(tjsch->max_cps_ram,&rsnaps);CHKERRQ(ierr); 2202 revolve_create_offline(rfine,rsnaps); 2203 break; 2204 case TWO_LEVEL_TWO_REVOLVE: 2205 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. */ 2206 diskstack->stacksize = diskblocks; 2207 ierr = PetscRevolveIntCast(tjsch->stride,&rfine);CHKERRQ(ierr); 2208 ierr = PetscRevolveIntCast(tjsch->max_cps_ram,&rsnaps);CHKERRQ(ierr); 2209 revolve_create_offline(rfine,rsnaps); 2210 ierr = PetscRevolveIntCast((tjsch->total_steps+tjsch->stride-1)/tjsch->stride,&rfine);CHKERRQ(ierr); 2211 ierr = PetscRevolveIntCast(diskblocks,&rsnaps);CHKERRQ(ierr); 2212 revolve2_create_offline(rfine,rsnaps); 2213 ierr = PetscNew(&rctx2);CHKERRQ(ierr); 2214 rctx2->snaps_in = rsnaps; 2215 rctx2->reverseonestep = PETSC_FALSE; 2216 rctx2->check = 0; 2217 rctx2->oldcapo = 0; 2218 rctx2->capo = 0; 2219 rctx2->info = 2; 2220 rctx2->fine = rfine; 2221 tjsch->rctx2 = rctx2; 2222 diskstack->top = -1; 2223 ierr = PetscMalloc1(diskstack->stacksize*sizeof(PetscInt),&diskstack->container);CHKERRQ(ierr); 2224 break; 2225 case REVOLVE_OFFLINE: 2226 ierr = PetscRevolveIntCast(tjsch->total_steps,&rfine);CHKERRQ(ierr); 2227 ierr = PetscRevolveIntCast(tjsch->max_cps_ram,&rsnaps);CHKERRQ(ierr); 2228 revolve_create_offline(rfine,rsnaps); 2229 break; 2230 case REVOLVE_ONLINE: 2231 stack->stacksize = tjsch->max_cps_ram; 2232 ierr = PetscRevolveIntCast(tjsch->max_cps_ram,&rsnaps);CHKERRQ(ierr); 2233 revolve_create_online(rsnaps); 2234 break; 2235 case REVOLVE_MULTISTAGE: 2236 ierr = PetscRevolveIntCast(tjsch->total_steps,&rfine);CHKERRQ(ierr); 2237 ierr = PetscRevolveIntCast(tjsch->max_cps_ram,&rsnaps);CHKERRQ(ierr); 2238 ierr = PetscRevolveIntCast(tjsch->max_cps_ram+tjsch->max_cps_disk,&rsnaps2);CHKERRQ(ierr); 2239 revolve_create_multistage(rfine,rsnaps2,rsnaps); 2240 break; 2241 default: 2242 break; 2243 } 2244 ierr = PetscNew(&rctx);CHKERRQ(ierr); 2245 ierr = PetscRevolveIntCast(tjsch->max_cps_ram,&rsnaps);CHKERRQ(ierr); 2246 rctx->snaps_in = rsnaps; /* for theta methods snaps_in=2*max_cps_ram */ 2247 rctx->reverseonestep = PETSC_FALSE; 2248 rctx->check = 0; 2249 rctx->oldcapo = 0; 2250 rctx->capo = 0; 2251 rctx->info = 2; 2252 if (tjsch->stride > 1) { 2253 ierr = PetscRevolveIntCast(tjsch->stride,&rfine);CHKERRQ(ierr); 2254 } else { 2255 ierr = PetscRevolveIntCast(tjsch->total_steps,&rfine);CHKERRQ(ierr); 2256 } 2257 rctx->fine = rfine; 2258 tjsch->rctx = rctx; 2259 if (tjsch->stype == REVOLVE_ONLINE) rctx->fine = -1; 2260 #endif 2261 } else { 2262 if (tjsch->stype == TWO_LEVEL_NOREVOLVE) stack->stacksize = tjsch->stride-1; /* need tjsch->stride-1 at most */ 2263 if (tjsch->stype == NONE) { 2264 if (fixedtimestep) stack->stacksize = stack->solution_only ? tjsch->total_steps : tjsch->total_steps-1; 2265 else { /* adaptive time step */ 2266 /* if max_cps_ram is not specified, use maximal allowed number of steps for stack size */ 2267 if (tjsch->max_cps_ram == -1) stack->stacksize = ts->max_steps < PETSC_MAX_INT ? ts->max_steps : 10000; 2268 tjsch->total_steps = stack->solution_only ? stack->stacksize : stack->stacksize+1; /* will be updated as time integration advances */ 2269 } 2270 } 2271 } 2272 2273 if ((tjsch->stype >= TWO_LEVEL_NOREVOLVE && tjsch->stype < REVOLVE_OFFLINE) || tjsch->stype == REVOLVE_MULTISTAGE) { /* these types need to use disk */ 2274 ierr = TSTrajectorySetUp_Basic(tj,ts);CHKERRQ(ierr); 2275 } 2276 2277 stack->stacksize = PetscMax(stack->stacksize,1); 2278 tjsch->recompute = PETSC_FALSE; 2279 ierr = StackInit(stack,stack->stacksize,numY);CHKERRQ(ierr); 2280 PetscFunctionReturn(0); 2281 } 2282 2283 static PetscErrorCode TSTrajectoryReset_Memory(TSTrajectory tj) 2284 { 2285 #if defined (PETSC_HAVE_REVOLVE) || defined (PETSC_HAVE_CAMS) 2286 TJScheduler *tjsch = (TJScheduler*)tj->data; 2287 PetscErrorCode ierr; 2288 #endif 2289 2290 PetscFunctionBegin; 2291 #if defined(PETSC_HAVE_REVOLVE) 2292 if (tjsch->stype > TWO_LEVEL_NOREVOLVE) { 2293 revolve_reset(); 2294 if (tjsch->stype == TWO_LEVEL_TWO_REVOLVE) { 2295 revolve2_reset(); 2296 ierr = PetscFree(tjsch->diskstack.container);CHKERRQ(ierr); 2297 } 2298 } 2299 if (tjsch->stype > TWO_LEVEL_NOREVOLVE) { 2300 ierr = PetscFree(tjsch->rctx);CHKERRQ(ierr); 2301 ierr = PetscFree(tjsch->rctx2);CHKERRQ(ierr); 2302 } 2303 #endif 2304 #if defined(PETSC_HAVE_CAMS) 2305 if (tjsch->stype == CAMS_OFFLINE) { 2306 if (tjsch->stack.solution_only) offline_ca_destroy(); 2307 else offline_ca_destroy(); 2308 ierr = PetscFree(tjsch->actx);CHKERRQ(ierr); 2309 } 2310 #endif 2311 PetscFunctionReturn(0); 2312 } 2313 2314 static PetscErrorCode TSTrajectoryDestroy_Memory(TSTrajectory tj) 2315 { 2316 TJScheduler *tjsch = (TJScheduler*)tj->data; 2317 PetscErrorCode ierr; 2318 2319 PetscFunctionBegin; 2320 ierr = StackDestroy(&tjsch->stack);CHKERRQ(ierr); 2321 ierr = PetscViewerDestroy(&tjsch->viewer);CHKERRQ(ierr); 2322 ierr = PetscObjectComposeFunction((PetscObject)tj,"TSTrajectorySetMaxCpsRAM_C",NULL);CHKERRQ(ierr); 2323 ierr = PetscObjectComposeFunction((PetscObject)tj,"TSTrajectorySetMaxCpsDisk_C",NULL);CHKERRQ(ierr); 2324 ierr = PetscObjectComposeFunction((PetscObject)tj,"TSTrajectorySetMaxUnitsRAM_C",NULL);CHKERRQ(ierr); 2325 ierr = PetscObjectComposeFunction((PetscObject)tj,"TSTrajectorySetMaxUnitsDisk_C",NULL);CHKERRQ(ierr); 2326 ierr = PetscFree(tjsch);CHKERRQ(ierr); 2327 PetscFunctionReturn(0); 2328 } 2329 2330 /*MC 2331 TSTRAJECTORYMEMORY - Stores each solution of the ODE/ADE in memory 2332 2333 Level: intermediate 2334 2335 .seealso: TSTrajectoryCreate(), TS, TSTrajectorySetType() 2336 2337 M*/ 2338 PETSC_EXTERN PetscErrorCode TSTrajectoryCreate_Memory(TSTrajectory tj,TS ts) 2339 { 2340 TJScheduler *tjsch; 2341 PetscErrorCode ierr; 2342 2343 PetscFunctionBegin; 2344 tj->ops->set = TSTrajectorySet_Memory; 2345 tj->ops->get = TSTrajectoryGet_Memory; 2346 tj->ops->setup = TSTrajectorySetUp_Memory; 2347 tj->ops->setfromoptions = TSTrajectorySetFromOptions_Memory; 2348 tj->ops->reset = TSTrajectoryReset_Memory; 2349 tj->ops->destroy = TSTrajectoryDestroy_Memory; 2350 2351 ierr = PetscNew(&tjsch);CHKERRQ(ierr); 2352 tjsch->stype = NONE; 2353 tjsch->max_cps_ram = -1; /* -1 indicates that it is not set */ 2354 tjsch->max_cps_disk = -1; /* -1 indicates that it is not set */ 2355 tjsch->stride = 0; /* if not zero, two-level checkpointing will be used */ 2356 #if defined(PETSC_HAVE_REVOLVE) 2357 tjsch->use_online = PETSC_FALSE; 2358 #endif 2359 tjsch->save_stack = PETSC_TRUE; 2360 2361 tjsch->stack.solution_only = tj->solution_only; 2362 ierr = PetscViewerCreate(PetscObjectComm((PetscObject)tj),&tjsch->viewer);CHKERRQ(ierr); 2363 ierr = PetscViewerSetType(tjsch->viewer,PETSCVIEWERBINARY);CHKERRQ(ierr); 2364 ierr = PetscViewerPushFormat(tjsch->viewer,PETSC_VIEWER_NATIVE);CHKERRQ(ierr); 2365 ierr = PetscViewerFileSetMode(tjsch->viewer,FILE_MODE_WRITE);CHKERRQ(ierr); 2366 2367 ierr = PetscObjectComposeFunction((PetscObject)tj,"TSTrajectorySetMaxCpsRAM_C",TSTrajectorySetMaxCpsRAM_Memory);CHKERRQ(ierr); 2368 ierr = PetscObjectComposeFunction((PetscObject)tj,"TSTrajectorySetMaxCpsDisk_C",TSTrajectorySetMaxCpsDisk_Memory);CHKERRQ(ierr); 2369 ierr = PetscObjectComposeFunction((PetscObject)tj,"TSTrajectorySetMaxUnitsRAM_C",TSTrajectorySetMaxUnitsRAM_Memory);CHKERRQ(ierr); 2370 ierr = PetscObjectComposeFunction((PetscObject)tj,"TSTrajectorySetMaxUnitsDisk_C",TSTrajectorySetMaxUnitsDisk_Memory);CHKERRQ(ierr); 2371 tj->data = tjsch; 2372 PetscFunctionReturn(0); 2373 } 2374