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