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