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