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