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