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 #endif 6 7 PetscLogEvent TSTrajectory_DiskWrite, TSTrajectory_DiskRead; 8 9 typedef enum {NONE,TWO_LEVEL_NOREVOLVE,TWO_LEVEL_REVOLVE,TWO_LEVEL_TWO_REVOLVE,REVOLVE_OFFLINE,REVOLVE_ONLINE,REVOLVE_MULTISTAGE} SchedulerType; 10 11 typedef struct _StackElement { 12 PetscInt stepnum; 13 Vec X; 14 Vec *Y; 15 PetscReal time; 16 PetscReal timeprev; /* for no solution_only mode */ 17 PetscReal timenext; /* for solution_only mode */ 18 } *StackElement; 19 20 #if defined(PETSC_HAVE_REVOLVE) 21 typedef struct _RevolveCTX { 22 PetscBool reverseonestep; 23 PetscInt where; 24 PetscInt snaps_in; 25 PetscInt stepsleft; 26 PetscInt check; 27 PetscInt oldcapo; 28 PetscInt capo; 29 PetscInt fine; 30 PetscInt info; 31 } RevolveCTX; 32 #endif 33 34 typedef struct _Stack { 35 PetscInt stacksize; 36 PetscInt top; 37 StackElement *container; 38 PetscInt numY; 39 PetscBool solution_only; 40 PetscBool use_dram; 41 } Stack; 42 43 typedef struct _DiskStack { 44 PetscInt stacksize; 45 PetscInt top; 46 PetscInt *container; 47 } DiskStack; 48 49 typedef struct _TJScheduler { 50 SchedulerType stype; 51 #if defined(PETSC_HAVE_REVOLVE) 52 RevolveCTX *rctx,*rctx2; 53 PetscBool use_online; 54 PetscInt store_stride; 55 #endif 56 PetscBool recompute; 57 PetscBool skip_trajectory; 58 PetscBool save_stack; 59 PetscInt max_cps_ram; /* maximum checkpoints in RAM */ 60 PetscInt max_cps_disk; /* maximum checkpoints on disk */ 61 PetscInt stride; 62 PetscInt total_steps; /* total number of steps */ 63 Stack stack; 64 DiskStack diskstack; 65 } TJScheduler; 66 67 static PetscErrorCode TurnForwardWithStepsize(TS ts,PetscReal nextstepsize) 68 { 69 PetscReal stepsize; 70 PetscErrorCode ierr; 71 72 PetscFunctionBegin; 73 /* reverse the direction */ 74 ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr); 75 stepsize = nextstepsize; 76 ierr = TSSetTimeStep(ts,stepsize);CHKERRQ(ierr); 77 PetscFunctionReturn(0); 78 } 79 80 static PetscErrorCode TurnForward(TS ts) 81 { 82 PetscReal stepsize; 83 PetscErrorCode ierr; 84 85 PetscFunctionBegin; 86 /* reverse the direction */ 87 ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr); 88 ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr); 89 PetscFunctionReturn(0); 90 } 91 92 static PetscErrorCode TurnBackward(TS ts) 93 { 94 PetscReal stepsize; 95 PetscErrorCode ierr; 96 97 PetscFunctionBegin; 98 /* reverse the direction */ 99 stepsize = ts->ptime_prev-ts->ptime; 100 ierr = TSSetTimeStep(ts,stepsize);CHKERRQ(ierr); 101 PetscFunctionReturn(0); 102 } 103 104 static PetscErrorCode ElementCreate(TS ts,Stack *stack,StackElement *e) 105 { 106 Vec X; 107 Vec *Y; 108 PetscErrorCode ierr; 109 110 PetscFunctionBegin; 111 if (stack->use_dram) { 112 ierr = PetscMallocSetDRAM();CHKERRQ(ierr); 113 } 114 ierr = PetscCalloc1(1,e);CHKERRQ(ierr); 115 ierr = TSGetSolution(ts,&X);CHKERRQ(ierr); 116 ierr = VecDuplicate(X,&(*e)->X);CHKERRQ(ierr); 117 if (stack->numY > 0 && !stack->solution_only) { 118 ierr = TSGetStages(ts,&stack->numY,&Y);CHKERRQ(ierr); 119 ierr = VecDuplicateVecs(Y[0],stack->numY,&(*e)->Y);CHKERRQ(ierr); 120 } 121 if (stack->use_dram) { 122 ierr = PetscMallocResetDRAM();CHKERRQ(ierr); 123 } 124 PetscFunctionReturn(0); 125 } 126 127 static PetscErrorCode ElementSet(TS ts,Stack *stack,StackElement *e,PetscInt stepnum,PetscReal time,Vec X) 128 { 129 Vec *Y; 130 PetscInt i; 131 PetscReal timeprev; 132 PetscErrorCode ierr; 133 134 PetscFunctionBegin; 135 ierr = VecCopy(X,(*e)->X);CHKERRQ(ierr); 136 if (stack->numY > 0 && !stack->solution_only) { 137 ierr = TSGetStages(ts,&stack->numY,&Y);CHKERRQ(ierr); 138 for (i=0;i<stack->numY;i++) { 139 ierr = VecCopy(Y[i],(*e)->Y[i]);CHKERRQ(ierr); 140 } 141 } 142 (*e)->stepnum = stepnum; 143 (*e)->time = time; 144 /* for consistency */ 145 if (stepnum == 0) { 146 (*e)->timeprev = (*e)->time - ts->time_step; 147 } else { 148 ierr = TSGetPrevTime(ts,&timeprev);CHKERRQ(ierr); 149 (*e)->timeprev = timeprev; 150 } 151 PetscFunctionReturn(0); 152 } 153 154 static PetscErrorCode ElementDestroy(Stack *stack,StackElement e) 155 { 156 PetscErrorCode ierr; 157 158 PetscFunctionBegin; 159 if (stack->use_dram) { 160 ierr = PetscMallocSetDRAM();CHKERRQ(ierr); 161 } 162 ierr = VecDestroy(&e->X);CHKERRQ(ierr); 163 if (stack->numY > 0 && !stack->solution_only) { 164 ierr = VecDestroyVecs(stack->numY,&e->Y);CHKERRQ(ierr); 165 } 166 ierr = PetscFree(e);CHKERRQ(ierr); 167 if (stack->use_dram) { 168 ierr = PetscMallocResetDRAM();CHKERRQ(ierr); 169 } 170 PetscFunctionReturn(0); 171 } 172 173 static PetscErrorCode StackResize(Stack *stack,PetscInt newsize) 174 { 175 StackElement *newcontainer; 176 PetscInt i; 177 PetscErrorCode ierr; 178 179 PetscFunctionBegin; 180 ierr = PetscMalloc1(newsize*sizeof(StackElement),&newcontainer);CHKERRQ(ierr); 181 for (i=0;i<stack->stacksize;i++) { 182 newcontainer[i] = stack->container[i]; 183 } 184 ierr = PetscFree(stack->container);CHKERRQ(ierr); 185 stack->container = newcontainer; 186 stack->stacksize = newsize; 187 PetscFunctionReturn(0); 188 } 189 190 static PetscErrorCode StackPush(Stack *stack,StackElement e) 191 { 192 PetscFunctionBegin; 193 if (stack->top+1 >= stack->stacksize) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MEMC,"Maximum stack size (%D) exceeded",stack->stacksize); 194 stack->container[++stack->top] = e; 195 PetscFunctionReturn(0); 196 } 197 198 static PetscErrorCode StackPop(Stack *stack,StackElement *e) 199 { 200 PetscFunctionBegin; 201 if (stack->top == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEMC,"Empty stack"); 202 *e = stack->container[stack->top--]; 203 PetscFunctionReturn(0); 204 } 205 206 static PetscErrorCode StackTop(Stack *stack,StackElement *e) 207 { 208 PetscFunctionBegin; 209 *e = stack->container[stack->top]; 210 PetscFunctionReturn(0); 211 } 212 213 static PetscErrorCode StackCreate(Stack *stack,PetscInt size,PetscInt ny) 214 { 215 PetscErrorCode ierr; 216 217 PetscFunctionBegin; 218 stack->top = -1; 219 stack->numY = ny; 220 221 ierr = PetscMalloc1(size*sizeof(StackElement),&stack->container);CHKERRQ(ierr); 222 PetscFunctionReturn(0); 223 } 224 225 static PetscErrorCode StackDestroy(Stack *stack) 226 { 227 PetscInt i,n; 228 StackElement e; 229 PetscErrorCode ierr; 230 231 PetscFunctionBegin; 232 if (stack->top > -1) { 233 n = stack->top+1; /* number of elements in the stack */ 234 for (i=0;i<n;i++) { 235 ierr = StackPop(stack,&e);CHKERRQ(ierr); 236 ierr = ElementDestroy(stack,e);CHKERRQ(ierr); 237 } 238 } 239 ierr = PetscFree(stack->container);CHKERRQ(ierr); 240 PetscFunctionReturn(0); 241 } 242 243 #if defined(PETSC_HAVE_REVOLVE) 244 static PetscErrorCode StackFind(Stack *stack,StackElement *e,PetscInt index) 245 { 246 PetscFunctionBegin; 247 *e = stack->container[index]; 248 PetscFunctionReturn(0); 249 } 250 #endif 251 252 static PetscErrorCode OutputBIN(MPI_Comm comm,const char *filename,PetscViewer *viewer) 253 { 254 PetscErrorCode ierr; 255 256 PetscFunctionBegin; 257 ierr = PetscViewerCreate(comm,viewer);CHKERRQ(ierr); 258 ierr = PetscViewerSetType(*viewer,PETSCVIEWERBINARY);CHKERRQ(ierr); 259 ierr = PetscViewerFileSetMode(*viewer,FILE_MODE_WRITE);CHKERRQ(ierr); 260 ierr = PetscViewerFileSetName(*viewer,filename);CHKERRQ(ierr); 261 PetscFunctionReturn(0); 262 } 263 264 static PetscErrorCode WriteToDisk(PetscInt stepnum,PetscReal time,PetscReal timeprev,Vec X,Vec *Y,PetscInt numY,PetscBool solution_only,PetscViewer viewer) 265 { 266 PetscInt i; 267 PetscErrorCode ierr; 268 269 PetscFunctionBegin; 270 ierr = PetscViewerBinaryWrite(viewer,&stepnum,1,PETSC_INT,PETSC_FALSE);CHKERRQ(ierr); 271 ierr = VecView(X,viewer);CHKERRQ(ierr); 272 for (i=0;!solution_only && i<numY;i++) { 273 ierr = VecView(Y[i],viewer);CHKERRQ(ierr); 274 } 275 ierr = PetscViewerBinaryWrite(viewer,&time,1,PETSC_REAL,PETSC_FALSE);CHKERRQ(ierr); 276 ierr = PetscViewerBinaryWrite(viewer,&timeprev,1,PETSC_REAL,PETSC_FALSE);CHKERRQ(ierr); 277 PetscFunctionReturn(0); 278 } 279 280 static PetscErrorCode ReadFromDisk(PetscInt *stepnum,PetscReal *time,PetscReal *timeprev,Vec X,Vec *Y,PetscInt numY,PetscBool solution_only,PetscViewer viewer) 281 { 282 PetscInt i; 283 PetscErrorCode ierr; 284 285 PetscFunctionBegin; 286 ierr = PetscViewerBinaryRead(viewer,stepnum,1,NULL,PETSC_INT);CHKERRQ(ierr); 287 ierr = VecLoad(X,viewer);CHKERRQ(ierr); 288 for (i=0;!solution_only && i<numY;i++) { 289 ierr = VecLoad(Y[i],viewer);CHKERRQ(ierr); 290 } 291 ierr = PetscViewerBinaryRead(viewer,time,1,NULL,PETSC_REAL);CHKERRQ(ierr); 292 ierr = PetscViewerBinaryRead(viewer,timeprev,1,NULL,PETSC_REAL);CHKERRQ(ierr); 293 PetscFunctionReturn(0); 294 } 295 296 static PetscErrorCode StackDumpAll(TSTrajectory tj,TS ts,Stack *stack,PetscInt id) 297 { 298 Vec *Y; 299 PetscInt i; 300 StackElement e = NULL; 301 PetscViewer viewer; 302 char filename[PETSC_MAX_PATH_LEN]; 303 PetscErrorCode ierr; 304 MPI_Comm comm; 305 306 PetscFunctionBegin; 307 ierr = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr); 308 if (tj->monitor) { 309 ierr = PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 310 ierr = PetscViewerASCIIPrintf(tj->monitor,"Dump stack to file\n");CHKERRQ(ierr); 311 ierr = PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 312 } 313 if (id == 1) { 314 PetscMPIInt rank; 315 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 316 if (!rank) { 317 ierr = PetscRMTree("SA-data");CHKERRQ(ierr); 318 ierr = PetscMkdir("SA-data");CHKERRQ(ierr); 319 } 320 } 321 ierr = PetscSNPrintf(filename,sizeof(filename),"SA-data/SA-STACK%06d.bin",id);CHKERRQ(ierr); 322 ierr = OutputBIN(comm,filename,&viewer);CHKERRQ(ierr); 323 for (i=0;i<stack->stacksize;i++) { 324 e = stack->container[i]; 325 ierr = PetscLogEventBegin(TSTrajectory_DiskWrite,ts,0,0,0);CHKERRQ(ierr); 326 ierr = WriteToDisk(e->stepnum,e->time,e->timeprev,e->X,e->Y,stack->numY,stack->solution_only,viewer);CHKERRQ(ierr); 327 ierr = PetscLogEventEnd(TSTrajectory_DiskWrite,ts,0,0,0);CHKERRQ(ierr); 328 ts->trajectory->diskwrites++; 329 } 330 /* 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 */ 331 ierr = TSGetStages(ts,&stack->numY,&Y);CHKERRQ(ierr); 332 ierr = PetscLogEventBegin(TSTrajectory_DiskWrite,ts,0,0,0);CHKERRQ(ierr); 333 ierr = WriteToDisk(ts->total_steps,ts->ptime,ts->ptime_prev,ts->vec_sol,Y,stack->numY,stack->solution_only,viewer);CHKERRQ(ierr); 334 ierr = PetscLogEventEnd(TSTrajectory_DiskWrite,ts,0,0,0);CHKERRQ(ierr); 335 ts->trajectory->diskwrites++; 336 for (i=0;i<stack->stacksize;i++) { 337 ierr = StackPop(stack,&e);CHKERRQ(ierr); 338 ierr = ElementDestroy(stack,e);CHKERRQ(ierr); 339 } 340 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 341 PetscFunctionReturn(0); 342 } 343 344 static PetscErrorCode StackLoadAll(TSTrajectory tj,TS ts,Stack *stack,PetscInt id) 345 { 346 Vec *Y; 347 PetscInt i; 348 StackElement e; 349 PetscViewer viewer; 350 char filename[PETSC_MAX_PATH_LEN]; 351 PetscErrorCode ierr; 352 353 PetscFunctionBegin; 354 if (tj->monitor) { 355 ierr = PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 356 ierr = PetscViewerASCIIPrintf(tj->monitor,"Load stack from file\n");CHKERRQ(ierr); 357 ierr = PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 358 } 359 ierr = PetscSNPrintf(filename,sizeof filename,"SA-data/SA-STACK%06d.bin",id);CHKERRQ(ierr); 360 ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);CHKERRQ(ierr); 361 for (i=0;i<stack->stacksize;i++) { 362 ierr = ElementCreate(ts,stack,&e);CHKERRQ(ierr); 363 ierr = StackPush(stack,e);CHKERRQ(ierr); 364 ierr = PetscLogEventBegin(TSTrajectory_DiskRead,ts,0,0,0);CHKERRQ(ierr); 365 ierr = ReadFromDisk(&e->stepnum,&e->time,&e->timeprev,e->X,e->Y,stack->numY,stack->solution_only,viewer);CHKERRQ(ierr); 366 ierr = PetscLogEventEnd(TSTrajectory_DiskRead,ts,0,0,0);CHKERRQ(ierr); 367 ts->trajectory->diskreads++; 368 } 369 /* load the last step into TS */ 370 ierr = TSGetStages(ts,&stack->numY,&Y);CHKERRQ(ierr); 371 ierr = PetscLogEventBegin(TSTrajectory_DiskRead,ts,0,0,0);CHKERRQ(ierr); 372 ierr = ReadFromDisk(&ts->total_steps,&ts->ptime,&ts->ptime_prev,ts->vec_sol,Y,stack->numY,stack->solution_only,viewer);CHKERRQ(ierr); 373 ierr = PetscLogEventEnd(TSTrajectory_DiskRead,ts,0,0,0);CHKERRQ(ierr); 374 ts->trajectory->diskreads++; 375 ierr = TurnBackward(ts);CHKERRQ(ierr); 376 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 377 PetscFunctionReturn(0); 378 } 379 380 #if defined(PETSC_HAVE_REVOLVE) 381 static PetscErrorCode StackLoadLast(TSTrajectory tj,TS ts,Stack *stack,PetscInt id) 382 { 383 Vec *Y; 384 PetscInt size; 385 PetscViewer viewer; 386 char filename[PETSC_MAX_PATH_LEN]; 387 #if defined(PETSC_HAVE_MPIIO) 388 PetscBool usempiio; 389 #endif 390 int fd; 391 off_t off,offset; 392 PetscErrorCode ierr; 393 394 PetscFunctionBegin; 395 if (tj->monitor) { 396 ierr = PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 397 ierr = PetscViewerASCIIPrintf(tj->monitor,"Load last stack element from file\n");CHKERRQ(ierr); 398 ierr = PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 399 } 400 ierr = TSGetStages(ts,&stack->numY,&Y);CHKERRQ(ierr); 401 ierr = VecGetSize(Y[0],&size);CHKERRQ(ierr); 402 /* VecView writes to file two extra int's for class id and number of rows */ 403 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; 404 405 ierr = PetscSNPrintf(filename,sizeof filename,"SA-data/SA-STACK%06d.bin",id);CHKERRQ(ierr); 406 ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);CHKERRQ(ierr); 407 #if defined(PETSC_HAVE_MPIIO) 408 ierr = PetscViewerBinaryGetUseMPIIO(viewer,&usempiio); 409 if (usempiio) { 410 ierr = PetscViewerBinaryGetMPIIODescriptor(viewer,(MPI_File*)&fd);CHKERRQ(ierr); 411 ierr = PetscBinarySynchronizedSeek(PETSC_COMM_WORLD,fd,off,PETSC_BINARY_SEEK_END,&offset);CHKERRQ(ierr); 412 } else { 413 #endif 414 ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 415 ierr = PetscBinarySeek(fd,off,PETSC_BINARY_SEEK_END,&offset);CHKERRQ(ierr); 416 #if defined(PETSC_HAVE_MPIIO) 417 } 418 #endif 419 /* load the last step into TS */ 420 ierr = PetscLogEventBegin(TSTrajectory_DiskRead,ts,0,0,0);CHKERRQ(ierr); 421 ierr = ReadFromDisk(&ts->total_steps,&ts->ptime,&ts->ptime_prev,ts->vec_sol,Y,stack->numY,stack->solution_only,viewer);CHKERRQ(ierr); 422 ierr = PetscLogEventEnd(TSTrajectory_DiskRead,ts,0,0,0);CHKERRQ(ierr); 423 ts->trajectory->diskreads++; 424 ierr = TurnBackward(ts);CHKERRQ(ierr); 425 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 426 PetscFunctionReturn(0); 427 } 428 #endif 429 430 static PetscErrorCode DumpSingle(TSTrajectory tj,TS ts,Stack *stack,PetscInt id) 431 { 432 Vec *Y; 433 PetscInt stepnum; 434 PetscViewer viewer; 435 char filename[PETSC_MAX_PATH_LEN]; 436 PetscErrorCode ierr; 437 MPI_Comm comm; 438 439 PetscFunctionBegin; 440 ierr = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr); 441 if (tj->monitor) { 442 ierr = PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 443 ierr = PetscViewerASCIIPrintf(tj->monitor,"Load a single point from file\n");CHKERRQ(ierr); 444 ierr = PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 445 } 446 ierr = TSGetTotalSteps(ts,&stepnum);CHKERRQ(ierr); 447 if (id == 1) { 448 PetscMPIInt rank; 449 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 450 if (!rank) { 451 ierr = PetscRMTree("SA-data");CHKERRQ(ierr); 452 ierr = PetscMkdir("SA-data");CHKERRQ(ierr); 453 } 454 } 455 ierr = PetscSNPrintf(filename,sizeof(filename),"SA-data/SA-CPS%06d.bin",id);CHKERRQ(ierr); 456 ierr = OutputBIN(comm,filename,&viewer);CHKERRQ(ierr); 457 458 ierr = TSGetStages(ts,&stack->numY,&Y);CHKERRQ(ierr); 459 ierr = PetscLogEventBegin(TSTrajectory_DiskWrite,ts,0,0,0);CHKERRQ(ierr); 460 ierr = WriteToDisk(stepnum,ts->ptime,ts->ptime_prev,ts->vec_sol,Y,stack->numY,stack->solution_only,viewer);CHKERRQ(ierr); 461 ierr = PetscLogEventEnd(TSTrajectory_DiskWrite,ts,0,0,0);CHKERRQ(ierr); 462 ts->trajectory->diskwrites++; 463 464 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 465 PetscFunctionReturn(0); 466 } 467 468 static PetscErrorCode LoadSingle(TSTrajectory tj,TS ts,Stack *stack,PetscInt id) 469 { 470 Vec *Y; 471 PetscViewer viewer; 472 char filename[PETSC_MAX_PATH_LEN]; 473 PetscErrorCode ierr; 474 475 PetscFunctionBegin; 476 if (tj->monitor) { 477 ierr = PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 478 ierr = PetscViewerASCIIPrintf(tj->monitor,"Load a single point from file\n");CHKERRQ(ierr); 479 ierr = PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 480 } 481 ierr = PetscSNPrintf(filename,sizeof filename,"SA-data/SA-CPS%06d.bin",id);CHKERRQ(ierr); 482 ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);CHKERRQ(ierr); 483 484 ierr = TSGetStages(ts,&stack->numY,&Y);CHKERRQ(ierr); 485 ierr = PetscLogEventBegin(TSTrajectory_DiskRead,ts,0,0,0);CHKERRQ(ierr); 486 ierr = ReadFromDisk(&ts->total_steps,&ts->ptime,&ts->ptime_prev,ts->vec_sol,Y,stack->numY,stack->solution_only,viewer);CHKERRQ(ierr); 487 ierr = PetscLogEventEnd(TSTrajectory_DiskRead,ts,0,0,0);CHKERRQ(ierr); 488 ts->trajectory->diskreads++; 489 490 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 491 PetscFunctionReturn(0); 492 } 493 494 static PetscErrorCode UpdateTS(TS ts,Stack *stack,StackElement e) 495 { 496 Vec *Y; 497 PetscInt i; 498 PetscErrorCode ierr; 499 500 PetscFunctionBegin; 501 ierr = VecCopy(e->X,ts->vec_sol);CHKERRQ(ierr); 502 if (!stack->solution_only) { 503 ierr = TSGetStages(ts,&stack->numY,&Y);CHKERRQ(ierr); 504 for (i=0;i<stack->numY;i++) { 505 ierr = VecCopy(e->Y[i],Y[i]);CHKERRQ(ierr); 506 } 507 } 508 ierr = TSSetTimeStep(ts,e->timeprev-e->time);CHKERRQ(ierr); /* stepsize will be negative */ 509 ts->ptime = e->time; 510 ts->ptime_prev = e->timeprev; 511 PetscFunctionReturn(0); 512 } 513 514 static PetscErrorCode ReCompute(TS ts,TJScheduler *tjsch,PetscInt stepnumbegin,PetscInt stepnumend) 515 { 516 Stack *stack = &tjsch->stack; 517 PetscInt i,adjsteps; 518 PetscErrorCode ierr; 519 520 PetscFunctionBegin; 521 adjsteps = ts->steps; 522 ts->steps = stepnumbegin; /* global step number */ 523 for (i=stepnumbegin;i<stepnumend;i++) { /* assume fixed step size */ 524 if (stack->solution_only && !tjsch->skip_trajectory) { /* revolve online need this */ 525 ierr = TSTrajectorySet(ts->trajectory,ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 526 } 527 ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 528 ierr = TSStep(ts);CHKERRQ(ierr); 529 if (!stack->solution_only && !tjsch->skip_trajectory) { 530 ierr = TSTrajectorySet(ts->trajectory,ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 531 } 532 ierr = TSEventHandler(ts);CHKERRQ(ierr); 533 if (!ts->steprollback) { 534 ierr = TSPostStep(ts);CHKERRQ(ierr); 535 } 536 } 537 ierr = TurnBackward(ts);CHKERRQ(ierr); 538 ts->trajectory->recomps += stepnumend-stepnumbegin; /* recomputation counter */ 539 ts->steps = adjsteps; 540 ts->total_steps = stepnumend; 541 PetscFunctionReturn(0); 542 } 543 544 static PetscErrorCode TopLevelStore(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscInt localstepnum,PetscInt laststridesize,PetscBool *done) 545 { 546 Stack *stack = &tjsch->stack; 547 DiskStack *diskstack = &tjsch->diskstack; 548 PetscInt stridenum; 549 PetscErrorCode ierr; 550 551 PetscFunctionBegin; 552 *done = PETSC_FALSE; 553 stridenum = stepnum/tjsch->stride; 554 /* make sure saved checkpoint id starts from 1 555 skip last stride when using stridenum+1 556 skip first stride when using stridenum */ 557 if (stack->solution_only) { 558 if (tjsch->save_stack) { 559 if (localstepnum == tjsch->stride-1 && stepnum < tjsch->total_steps-laststridesize) { /* current step will be saved without going through stack */ 560 ierr = StackDumpAll(tj,ts,stack,stridenum+1);CHKERRQ(ierr); 561 if (tjsch->stype == TWO_LEVEL_TWO_REVOLVE) diskstack->container[++diskstack->top] = stridenum+1; 562 *done = PETSC_TRUE; 563 } 564 } else { 565 if (localstepnum == 0 && stepnum < tjsch->total_steps-laststridesize) { 566 ierr = DumpSingle(tj,ts,stack,stridenum+1);CHKERRQ(ierr); 567 if (tjsch->stype == TWO_LEVEL_TWO_REVOLVE) diskstack->container[++diskstack->top] = stridenum+1; 568 *done = PETSC_TRUE; 569 } 570 } 571 } else { 572 if (tjsch->save_stack) { 573 if (localstepnum == 0 && stepnum < tjsch->total_steps && stepnum != 0) { /* skip the first stride */ 574 ierr = StackDumpAll(tj,ts,stack,stridenum);CHKERRQ(ierr); 575 if (tjsch->stype == TWO_LEVEL_TWO_REVOLVE) diskstack->container[++diskstack->top] = stridenum; 576 *done = PETSC_TRUE; 577 } 578 } else { 579 if (localstepnum == 1 && stepnum < tjsch->total_steps-laststridesize) { 580 ierr = DumpSingle(tj,ts,stack,stridenum+1);CHKERRQ(ierr); 581 if (tjsch->stype == TWO_LEVEL_TWO_REVOLVE) diskstack->container[++diskstack->top] = stridenum+1; 582 *done = PETSC_TRUE; 583 } 584 } 585 } 586 PetscFunctionReturn(0); 587 } 588 589 static PetscErrorCode SetTrajN(TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X) 590 { 591 Stack *stack = &tjsch->stack; 592 StackElement e; 593 PetscErrorCode ierr; 594 595 PetscFunctionBegin; 596 /* skip the last step */ 597 if (ts->reason) { /* only affect the forward run */ 598 /* update total_steps in the end of forward run */ 599 if (stepnum != tjsch->total_steps) tjsch->total_steps = stepnum; 600 if (stack->solution_only) { 601 /* get rid of the solution at second last step */ 602 ierr = StackPop(stack,&e);CHKERRQ(ierr); 603 ierr = ElementDestroy(stack,e);CHKERRQ(ierr); 604 } 605 PetscFunctionReturn(0); 606 } 607 /* do not save trajectory at the recompute stage for solution_only mode */ 608 if (stack->solution_only && tjsch->recompute) PetscFunctionReturn(0); 609 /* skip the first step for no_solution_only mode */ 610 if (!stack->solution_only && stepnum == 0) PetscFunctionReturn(0); 611 612 /* resize the stack */ 613 if (stack->top+1 == stack->stacksize) { 614 ierr = StackResize(stack,2*stack->stacksize);CHKERRQ(ierr); 615 } 616 /* update timenext for the previous step; necessary for step adaptivity */ 617 if (stack->top > -1) { 618 ierr = StackTop(stack,&e);CHKERRQ(ierr); 619 e->timenext = ts->ptime; 620 } 621 if (stepnum < stack->top) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_MEMC,"Illegal modification of a non-top stack element"); 622 ierr = ElementCreate(ts,stack,&e);CHKERRQ(ierr); 623 ierr = ElementSet(ts,stack,&e,stepnum,time,X);CHKERRQ(ierr); 624 ierr = StackPush(stack,e);CHKERRQ(ierr); 625 PetscFunctionReturn(0); 626 } 627 628 static PetscErrorCode GetTrajN(TS ts,TJScheduler *tjsch,PetscInt stepnum) 629 { 630 Stack *stack = &tjsch->stack; 631 StackElement e; 632 PetscErrorCode ierr; 633 634 PetscFunctionBegin; 635 if (stepnum == tjsch->total_steps) { 636 ierr = TurnBackward(ts);CHKERRQ(ierr); 637 PetscFunctionReturn(0); 638 } 639 /* restore a checkpoint */ 640 ierr = StackTop(stack,&e);CHKERRQ(ierr); 641 ierr = UpdateTS(ts,stack,e);CHKERRQ(ierr); 642 if (stack->solution_only) {/* recompute one step */ 643 tjsch->recompute = PETSC_TRUE; 644 ierr = TurnForwardWithStepsize(ts,e->timenext-e->time);CHKERRQ(ierr); 645 ierr = ReCompute(ts,tjsch,e->stepnum,stepnum);CHKERRQ(ierr); 646 } 647 ierr = StackPop(stack,&e);CHKERRQ(ierr); 648 ierr = ElementDestroy(stack,e);CHKERRQ(ierr); 649 PetscFunctionReturn(0); 650 } 651 652 static PetscErrorCode SetTrajTLNR(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X) 653 { 654 Stack *stack = &tjsch->stack; 655 PetscInt localstepnum,laststridesize; 656 StackElement e; 657 PetscBool done; 658 PetscErrorCode ierr; 659 660 PetscFunctionBegin; 661 if (!stack->solution_only && stepnum == 0) PetscFunctionReturn(0); 662 if (stack->solution_only && stepnum == tjsch->total_steps) PetscFunctionReturn(0); 663 if (tjsch->save_stack && tjsch->recompute) PetscFunctionReturn(0); 664 665 localstepnum = stepnum%tjsch->stride; 666 /* (stride size-1) checkpoints are saved in each stride; an extra point is added by StackDumpAll() */ 667 laststridesize = tjsch->total_steps%tjsch->stride; 668 if (!laststridesize) laststridesize = tjsch->stride; 669 670 if (!tjsch->recompute) { 671 ierr = TopLevelStore(tj,ts,tjsch,stepnum,localstepnum,laststridesize,&done);CHKERRQ(ierr); 672 if (!tjsch->save_stack && stepnum < tjsch->total_steps-laststridesize) PetscFunctionReturn(0); 673 } 674 if (!stack->solution_only && localstepnum == 0) PetscFunctionReturn(0); /* skip last point in each stride at recompute stage or last stride */ 675 if (stack->solution_only && localstepnum == tjsch->stride-1) PetscFunctionReturn(0); /* skip last step in each stride at recompute stage or last stride */ 676 677 ierr = ElementCreate(ts,stack,&e);CHKERRQ(ierr); 678 ierr = ElementSet(ts,stack,&e,stepnum,time,X);CHKERRQ(ierr); 679 ierr = StackPush(stack,e);CHKERRQ(ierr); 680 PetscFunctionReturn(0); 681 } 682 683 static PetscErrorCode GetTrajTLNR(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum) 684 { 685 Stack *stack = &tjsch->stack; 686 PetscInt id,localstepnum,laststridesize; 687 StackElement e; 688 PetscErrorCode ierr; 689 690 PetscFunctionBegin; 691 if (stepnum == tjsch->total_steps) { 692 ierr = TurnBackward(ts);CHKERRQ(ierr); 693 PetscFunctionReturn(0); 694 } 695 696 localstepnum = stepnum%tjsch->stride; 697 laststridesize = tjsch->total_steps%tjsch->stride; 698 if (!laststridesize) laststridesize = tjsch->stride; 699 if (stack->solution_only) { 700 /* fill stack with info */ 701 if (localstepnum == 0 && tjsch->total_steps-stepnum >= laststridesize) { 702 id = stepnum/tjsch->stride; 703 if (tjsch->save_stack) { 704 ierr = StackLoadAll(tj,ts,stack,id);CHKERRQ(ierr); 705 tjsch->recompute = PETSC_TRUE; 706 tjsch->skip_trajectory = PETSC_TRUE; 707 ierr = TurnForward(ts);CHKERRQ(ierr); 708 ierr = ReCompute(ts,tjsch,id*tjsch->stride-1,id*tjsch->stride);CHKERRQ(ierr); 709 tjsch->skip_trajectory = PETSC_FALSE; 710 } else { 711 ierr = LoadSingle(tj,ts,stack,id);CHKERRQ(ierr); 712 tjsch->recompute = PETSC_TRUE; 713 ierr = TurnForward(ts);CHKERRQ(ierr); 714 ierr = ReCompute(ts,tjsch,(id-1)*tjsch->stride,id*tjsch->stride);CHKERRQ(ierr); 715 } 716 PetscFunctionReturn(0); 717 } 718 /* restore a checkpoint */ 719 ierr = StackPop(stack,&e);CHKERRQ(ierr); 720 ierr = UpdateTS(ts,stack,e);CHKERRQ(ierr); 721 tjsch->recompute = PETSC_TRUE; 722 tjsch->skip_trajectory = PETSC_TRUE; 723 ierr = TurnForward(ts);CHKERRQ(ierr); 724 ierr = ReCompute(ts,tjsch,e->stepnum,stepnum);CHKERRQ(ierr); 725 tjsch->skip_trajectory = PETSC_FALSE; 726 ierr = ElementDestroy(stack,e);CHKERRQ(ierr); 727 } else { 728 /* fill stack with info */ 729 if (localstepnum == 0 && tjsch->total_steps-stepnum >= laststridesize) { 730 id = stepnum/tjsch->stride; 731 if (tjsch->save_stack) { 732 ierr = StackLoadAll(tj,ts,stack,id);CHKERRQ(ierr); 733 } else { 734 ierr = LoadSingle(tj,ts,stack,id);CHKERRQ(ierr); 735 ierr = ElementCreate(ts,stack,&e);CHKERRQ(ierr); 736 ierr = ElementSet(ts,stack,&e,(id-1)*tjsch->stride+1,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 737 ierr = StackPush(stack,e);CHKERRQ(ierr); 738 tjsch->recompute = PETSC_TRUE; 739 ierr = TurnForward(ts);CHKERRQ(ierr); 740 ierr = ReCompute(ts,tjsch,e->stepnum,id*tjsch->stride);CHKERRQ(ierr); 741 } 742 PetscFunctionReturn(0); 743 } 744 /* restore a checkpoint */ 745 ierr = StackPop(stack,&e);CHKERRQ(ierr); 746 ierr = UpdateTS(ts,stack,e);CHKERRQ(ierr); 747 ierr = ElementDestroy(stack,e);CHKERRQ(ierr); 748 } 749 PetscFunctionReturn(0); 750 } 751 752 #if defined(PETSC_HAVE_REVOLVE) 753 static PetscErrorCode printwhattodo(PetscViewer viewer,PetscInt whattodo,RevolveCTX *rctx,PetscInt shift) 754 { 755 PetscErrorCode ierr; 756 757 PetscFunctionBegin; 758 if (!viewer) PetscFunctionReturn(0); 759 760 switch(whattodo) { 761 case 1: 762 ierr = PetscViewerASCIIPrintf(viewer,"Advance from %D to %D\n",rctx->oldcapo+shift,rctx->capo+shift);CHKERRQ(ierr); 763 break; 764 case 2: 765 ierr = PetscViewerASCIIPrintf(viewer,"Store in checkpoint number %D (located in RAM)\n",rctx->check);CHKERRQ(ierr); 766 break; 767 case 3: 768 ierr = PetscViewerASCIIPrintf(viewer,"First turn: Initialize adjoints and reverse first step\n");CHKERRQ(ierr); 769 break; 770 case 4: 771 ierr = PetscViewerASCIIPrintf(viewer,"Forward and reverse one step\n");CHKERRQ(ierr); 772 break; 773 case 5: 774 ierr = PetscViewerASCIIPrintf(viewer,"Restore in checkpoint number %D (located in RAM)\n",rctx->check);CHKERRQ(ierr); 775 break; 776 case 7: 777 ierr = PetscViewerASCIIPrintf(viewer,"Store in checkpoint number %D (located on disk)\n",rctx->check);CHKERRQ(ierr); 778 break; 779 case 8: 780 ierr = PetscViewerASCIIPrintf(viewer,"Restore in checkpoint number %D (located on disk)\n",rctx->check);CHKERRQ(ierr); 781 break; 782 case -1: 783 ierr = PetscViewerASCIIPrintf(viewer,"Error!");CHKERRQ(ierr); 784 break; 785 } 786 PetscFunctionReturn(0); 787 } 788 789 static PetscErrorCode printwhattodo2(PetscViewer viewer,PetscInt whattodo,RevolveCTX *rctx,PetscInt shift) 790 { 791 PetscErrorCode ierr; 792 793 PetscFunctionBegin; 794 if (!viewer) PetscFunctionReturn(0); 795 796 switch(whattodo) { 797 case 1: 798 ierr = PetscViewerASCIIPrintf(viewer,"[Top Level] Advance from stride %D to stride %D\n",rctx->oldcapo+shift,rctx->capo+shift);CHKERRQ(ierr); 799 break; 800 case 2: 801 ierr = PetscViewerASCIIPrintf(viewer,"[Top Level] Store in checkpoint number %D\n",rctx->check);CHKERRQ(ierr); 802 break; 803 case 3: 804 ierr = PetscViewerASCIIPrintf(viewer,"[Top Level] First turn: Initialize adjoints and reverse first stride\n");CHKERRQ(ierr); 805 break; 806 case 4: 807 ierr = PetscViewerASCIIPrintf(viewer,"[Top Level] Forward and reverse one stride\n");CHKERRQ(ierr); 808 break; 809 case 5: 810 ierr = PetscViewerASCIIPrintf(viewer,"[Top Level] Restore in checkpoint number %D\n",rctx->check);CHKERRQ(ierr); 811 break; 812 case 7: 813 ierr = PetscViewerASCIIPrintf(viewer,"[Top Level] Store in top-level checkpoint number %D\n",rctx->check);CHKERRQ(ierr); 814 break; 815 case 8: 816 ierr = PetscViewerASCIIPrintf(viewer,"[Top Level] Restore in top-level checkpoint number %D\n",rctx->check);CHKERRQ(ierr); 817 break; 818 case -1: 819 ierr = PetscViewerASCIIPrintf(viewer,"[Top Level] Error!");CHKERRQ(ierr); 820 break; 821 } 822 PetscFunctionReturn(0); 823 } 824 825 static PetscErrorCode InitRevolve(PetscInt fine,PetscInt snaps,RevolveCTX *rctx) 826 { 827 PetscFunctionBegin; 828 revolve_reset(); 829 revolve_create_offline(fine,snaps); 830 rctx->snaps_in = snaps; 831 rctx->fine = fine; 832 rctx->check = 0; 833 rctx->capo = 0; 834 rctx->reverseonestep = PETSC_FALSE; 835 /* check stepsleft? */ 836 PetscFunctionReturn(0); 837 } 838 839 static PetscErrorCode FastForwardRevolve(RevolveCTX *rctx) 840 { 841 PetscInt whattodo; 842 843 PetscFunctionBegin; 844 whattodo = 0; 845 while(whattodo!=3) { /* we have to fast forward revolve to the beginning of the backward sweep due to unfriendly revolve interface */ 846 whattodo = revolve_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where); 847 } 848 PetscFunctionReturn(0); 849 } 850 851 static PetscErrorCode ApplyRevolve(PetscViewer viewer,SchedulerType stype,RevolveCTX *rctx,PetscInt total_steps,PetscInt stepnum,PetscInt localstepnum,PetscBool toplevel,PetscInt *store) 852 { 853 PetscErrorCode ierr; 854 PetscInt shift,whattodo; 855 856 PetscFunctionBegin; 857 *store = 0; 858 if (rctx->stepsleft > 0) { /* advance the solution without checkpointing anything as Revolve requires */ 859 rctx->stepsleft--; 860 PetscFunctionReturn(0); 861 } 862 /* let Revolve determine what to do next */ 863 shift = stepnum-localstepnum; 864 rctx->oldcapo = rctx->capo; 865 rctx->capo = localstepnum; 866 867 if (!toplevel) whattodo = revolve_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where); 868 else whattodo = revolve2_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where); 869 if (stype == REVOLVE_ONLINE && whattodo == 8) whattodo = 5; 870 if (stype == REVOLVE_ONLINE && whattodo == 7) whattodo = 2; 871 if (!toplevel) {ierr = printwhattodo(viewer,whattodo,rctx,shift);CHKERRQ(ierr);} 872 else {ierr = printwhattodo2(viewer,whattodo,rctx,shift);CHKERRQ(ierr);} 873 if (whattodo == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in the Revolve library"); 874 if (whattodo == 1) { /* advance some time steps */ 875 if (stype == REVOLVE_ONLINE && rctx->capo >= total_steps-1) { 876 revolve_turn(total_steps,&rctx->capo,&rctx->fine); 877 if (!toplevel) {ierr = printwhattodo(viewer,whattodo,rctx,shift);CHKERRQ(ierr);} 878 else {ierr = printwhattodo2(viewer,whattodo,rctx,shift);CHKERRQ(ierr);} 879 } 880 rctx->stepsleft = rctx->capo-rctx->oldcapo-1; 881 } 882 if (whattodo == 3 || whattodo == 4) { /* ready for a reverse step */ 883 rctx->reverseonestep = PETSC_TRUE; 884 } 885 if (whattodo == 5) { /* restore a checkpoint and ask Revolve what to do next */ 886 rctx->oldcapo = rctx->capo; 887 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*/ 888 else whattodo = revolve2_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where); 889 if (!toplevel) {ierr = printwhattodo(viewer,whattodo,rctx,shift);CHKERRQ(ierr);} 890 else {ierr = printwhattodo2(viewer,whattodo,rctx,shift);CHKERRQ(ierr);} 891 if (whattodo == 3 || whattodo == 4) rctx->reverseonestep = PETSC_TRUE; 892 if (whattodo == 1) rctx->stepsleft = rctx->capo-rctx->oldcapo; 893 } 894 if (whattodo == 7) { /* save the checkpoint to disk */ 895 *store = 2; 896 rctx->oldcapo = rctx->capo; 897 whattodo = revolve_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where); /* must return 1 */ 898 ierr = printwhattodo(viewer,whattodo,rctx,shift);CHKERRQ(ierr); 899 rctx->stepsleft = rctx->capo-rctx->oldcapo-1; 900 } 901 if (whattodo == 2) { /* store a checkpoint to RAM and ask Revolve how many time steps to advance next */ 902 *store = 1; 903 rctx->oldcapo = rctx->capo; 904 if (!toplevel) whattodo = revolve_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where); /* must return 1 */ 905 else whattodo = revolve2_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where); 906 if (!toplevel) {ierr = printwhattodo(viewer,whattodo,rctx,shift);CHKERRQ(ierr);} 907 else {ierr = printwhattodo2(viewer,whattodo,rctx,shift);CHKERRQ(ierr);} 908 if (stype == REVOLVE_ONLINE && rctx->capo >= total_steps-1) { 909 revolve_turn(total_steps,&rctx->capo,&rctx->fine); 910 ierr = printwhattodo(viewer,whattodo,rctx,shift);CHKERRQ(ierr); 911 } 912 rctx->stepsleft = rctx->capo-rctx->oldcapo-1; 913 } 914 PetscFunctionReturn(0); 915 } 916 917 static PetscErrorCode SetTrajROF(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X) 918 { 919 Stack *stack = &tjsch->stack; 920 PetscInt store; 921 StackElement e; 922 PetscErrorCode ierr; 923 924 PetscFunctionBegin; 925 if (!stack->solution_only && stepnum == 0) PetscFunctionReturn(0); 926 if (stack->solution_only && stepnum == tjsch->total_steps) PetscFunctionReturn(0); 927 ierr = ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,tjsch->total_steps,stepnum,stepnum,PETSC_FALSE,&store);CHKERRQ(ierr); 928 if (store == 1) { 929 if (stepnum < stack->top) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_MEMC,"Illegal modification of a non-top stack element"); 930 ierr = ElementCreate(ts,stack,&e);CHKERRQ(ierr); 931 ierr = ElementSet(ts,stack,&e,stepnum,time,X);CHKERRQ(ierr); 932 ierr = StackPush(stack,e);CHKERRQ(ierr); 933 } 934 PetscFunctionReturn(0); 935 } 936 937 static PetscErrorCode GetTrajROF(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum) 938 { 939 Stack *stack = &tjsch->stack; 940 PetscInt whattodo,shift,store; 941 StackElement e; 942 PetscErrorCode ierr; 943 944 PetscFunctionBegin; 945 if (stepnum == 0 || stepnum == tjsch->total_steps) { 946 ierr = TurnBackward(ts);CHKERRQ(ierr); 947 tjsch->rctx->reverseonestep = PETSC_FALSE; 948 PetscFunctionReturn(0); 949 } 950 /* restore a checkpoint */ 951 ierr = StackTop(stack,&e);CHKERRQ(ierr); 952 ierr = UpdateTS(ts,stack,e);CHKERRQ(ierr); 953 if (stack->solution_only) { /* start with restoring a checkpoint */ 954 tjsch->rctx->capo = stepnum; 955 tjsch->rctx->oldcapo = tjsch->rctx->capo; 956 shift = 0; 957 whattodo = revolve_action(&tjsch->rctx->check,&tjsch->rctx->capo,&tjsch->rctx->fine,tjsch->rctx->snaps_in,&tjsch->rctx->info,&tjsch->rctx->where); 958 ierr = printwhattodo(tj->monitor,whattodo,tjsch->rctx,shift);CHKERRQ(ierr); 959 } else { /* 2 revolve actions: restore a checkpoint and then advance */ 960 ierr = ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,tjsch->total_steps,stepnum,stepnum,PETSC_FALSE,&store);CHKERRQ(ierr); 961 if (tj->monitor) { 962 ierr = PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 963 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); 964 ierr = PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 965 } 966 if (!tjsch->rctx->reverseonestep && tjsch->rctx->stepsleft > 0) tjsch->rctx->stepsleft--; 967 } 968 if (stack->solution_only || (!stack->solution_only && e->stepnum < stepnum)) { 969 tjsch->recompute = PETSC_TRUE; 970 ierr = TurnForward(ts);CHKERRQ(ierr); 971 ierr = ReCompute(ts,tjsch,e->stepnum,stepnum);CHKERRQ(ierr); 972 } 973 if ((stack->solution_only && e->stepnum+1 == stepnum) || (!stack->solution_only && e->stepnum == stepnum)) { 974 ierr = StackPop(stack,&e);CHKERRQ(ierr); 975 ierr = ElementDestroy(stack,e);CHKERRQ(ierr); 976 } 977 tjsch->rctx->reverseonestep = PETSC_FALSE; 978 PetscFunctionReturn(0); 979 } 980 981 static PetscErrorCode SetTrajRON(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X) 982 { 983 Stack *stack = &tjsch->stack; 984 Vec *Y; 985 PetscInt i,store; 986 PetscReal timeprev; 987 StackElement e; 988 RevolveCTX *rctx = tjsch->rctx; 989 PetscErrorCode ierr; 990 991 PetscFunctionBegin; 992 if (!stack->solution_only && stepnum == 0) PetscFunctionReturn(0); 993 if (stack->solution_only && stepnum == tjsch->total_steps) PetscFunctionReturn(0); 994 ierr = ApplyRevolve(tj->monitor,tjsch->stype,rctx,tjsch->total_steps,stepnum,stepnum,PETSC_FALSE,&store);CHKERRQ(ierr); 995 if (store == 1) { 996 if (rctx->check != stack->top+1) { /* overwrite some non-top checkpoint in the stack */ 997 ierr = StackFind(stack,&e,rctx->check);CHKERRQ(ierr); 998 ierr = VecCopy(X,e->X);CHKERRQ(ierr); 999 if (stack->numY > 0 && !stack->solution_only) { 1000 ierr = TSGetStages(ts,&stack->numY,&Y);CHKERRQ(ierr); 1001 for (i=0;i<stack->numY;i++) { 1002 ierr = VecCopy(Y[i],e->Y[i]);CHKERRQ(ierr); 1003 } 1004 } 1005 e->stepnum = stepnum; 1006 e->time = time; 1007 ierr = TSGetPrevTime(ts,&timeprev);CHKERRQ(ierr); 1008 e->timeprev = timeprev; 1009 } else { 1010 if (stepnum < stack->top) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_MEMC,"Illegal modification of a non-top stack element"); 1011 ierr = ElementCreate(ts,stack,&e);CHKERRQ(ierr); 1012 ierr = ElementSet(ts,stack,&e,stepnum,time,X);CHKERRQ(ierr); 1013 ierr = StackPush(stack,e);CHKERRQ(ierr); 1014 } 1015 } 1016 PetscFunctionReturn(0); 1017 } 1018 1019 static PetscErrorCode GetTrajRON(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum) 1020 { 1021 Stack *stack = &tjsch->stack; 1022 PetscInt whattodo,shift; 1023 StackElement e; 1024 PetscErrorCode ierr; 1025 1026 PetscFunctionBegin; 1027 if (stepnum == 0 || stepnum == tjsch->total_steps) { 1028 ierr = TurnBackward(ts);CHKERRQ(ierr); 1029 tjsch->rctx->reverseonestep = PETSC_FALSE; 1030 PetscFunctionReturn(0); 1031 } 1032 tjsch->rctx->capo = stepnum; 1033 tjsch->rctx->oldcapo = tjsch->rctx->capo; 1034 shift = 0; 1035 whattodo = revolve_action(&tjsch->rctx->check,&tjsch->rctx->capo,&tjsch->rctx->fine,tjsch->rctx->snaps_in,&tjsch->rctx->info,&tjsch->rctx->where); /* whattodo=restore */ 1036 if (whattodo == 8) whattodo = 5; 1037 ierr = printwhattodo(tj->monitor,whattodo,tjsch->rctx,shift);CHKERRQ(ierr); 1038 /* restore a checkpoint */ 1039 ierr = StackFind(stack,&e,tjsch->rctx->check);CHKERRQ(ierr); 1040 ierr = UpdateTS(ts,stack,e);CHKERRQ(ierr); 1041 if (!stack->solution_only) { /* whattodo must be 5 */ 1042 /* ask Revolve what to do next */ 1043 tjsch->rctx->oldcapo = tjsch->rctx->capo; 1044 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*/ 1045 ierr = printwhattodo(tj->monitor,whattodo,tjsch->rctx,shift);CHKERRQ(ierr); 1046 if (whattodo == 3 || whattodo == 4) tjsch->rctx->reverseonestep = PETSC_TRUE; 1047 if (whattodo == 1) tjsch->rctx->stepsleft = tjsch->rctx->capo-tjsch->rctx->oldcapo; 1048 if (tj->monitor) { 1049 ierr = PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 1050 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); 1051 ierr = PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 1052 } 1053 if (!tjsch->rctx->reverseonestep && tjsch->rctx->stepsleft > 0) tjsch->rctx->stepsleft--; 1054 } 1055 if (stack->solution_only || (!stack->solution_only && e->stepnum < stepnum)) { 1056 tjsch->recompute = PETSC_TRUE; 1057 ierr = TurnForward(ts);CHKERRQ(ierr); 1058 ierr = ReCompute(ts,tjsch,e->stepnum,stepnum);CHKERRQ(ierr); 1059 } 1060 tjsch->rctx->reverseonestep = PETSC_FALSE; 1061 PetscFunctionReturn(0); 1062 } 1063 1064 static PetscErrorCode SetTrajTLR(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X) 1065 { 1066 Stack *stack = &tjsch->stack; 1067 PetscInt store,localstepnum,laststridesize; 1068 StackElement e; 1069 PetscBool done = PETSC_FALSE; 1070 PetscErrorCode ierr; 1071 1072 PetscFunctionBegin; 1073 if (!stack->solution_only && stepnum == 0) PetscFunctionReturn(0); 1074 if (stack->solution_only && stepnum == tjsch->total_steps) PetscFunctionReturn(0); 1075 1076 localstepnum = stepnum%tjsch->stride; 1077 laststridesize = tjsch->total_steps%tjsch->stride; 1078 if (!laststridesize) laststridesize = tjsch->stride; 1079 1080 if (!tjsch->recompute) { 1081 ierr = TopLevelStore(tj,ts,tjsch,stepnum,localstepnum,laststridesize,&done);CHKERRQ(ierr); 1082 /* revolve is needed for the last stride; different starting points for last stride between solutin_only and !solutin_only */ 1083 if (!stack->solution_only && !tjsch->save_stack && stepnum <= tjsch->total_steps-laststridesize) PetscFunctionReturn(0); 1084 if (stack->solution_only && !tjsch->save_stack && stepnum < tjsch->total_steps-laststridesize) PetscFunctionReturn(0); 1085 } 1086 if (tjsch->save_stack && done) { 1087 ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr); 1088 PetscFunctionReturn(0); 1089 } 1090 if (laststridesize < tjsch->stride) { 1091 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 */ 1092 ierr = InitRevolve(laststridesize,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr); 1093 } 1094 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 */ 1095 ierr = InitRevolve(laststridesize,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr); 1096 } 1097 } 1098 ierr = ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,tjsch->total_steps,stepnum,localstepnum,PETSC_FALSE,&store);CHKERRQ(ierr); 1099 if (store == 1) { 1100 if (localstepnum < stack->top) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_MEMC,"Illegal modification of a non-top stack element"); 1101 ierr = ElementCreate(ts,stack,&e);CHKERRQ(ierr); 1102 ierr = ElementSet(ts,stack,&e,stepnum,time,X);CHKERRQ(ierr); 1103 ierr = StackPush(stack,e);CHKERRQ(ierr); 1104 } 1105 PetscFunctionReturn(0); 1106 } 1107 1108 static PetscErrorCode GetTrajTLR(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum) 1109 { 1110 Stack *stack = &tjsch->stack; 1111 PetscInt whattodo,shift; 1112 PetscInt localstepnum,stridenum,laststridesize,store; 1113 StackElement e; 1114 PetscErrorCode ierr; 1115 1116 PetscFunctionBegin; 1117 localstepnum = stepnum%tjsch->stride; 1118 stridenum = stepnum/tjsch->stride; 1119 if (stepnum == tjsch->total_steps) { 1120 ierr = TurnBackward(ts);CHKERRQ(ierr); 1121 tjsch->rctx->reverseonestep = PETSC_FALSE; 1122 PetscFunctionReturn(0); 1123 } 1124 laststridesize = tjsch->total_steps%tjsch->stride; 1125 if (!laststridesize) laststridesize = tjsch->stride; 1126 if (stack->solution_only) { 1127 /* fill stack */ 1128 if (localstepnum == 0 && stepnum <= tjsch->total_steps-laststridesize) { 1129 if (tjsch->save_stack) { 1130 ierr = StackLoadAll(tj,ts,stack,stridenum);CHKERRQ(ierr); 1131 ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr); 1132 ierr = FastForwardRevolve(tjsch->rctx);CHKERRQ(ierr); 1133 tjsch->recompute = PETSC_TRUE; 1134 tjsch->skip_trajectory = PETSC_TRUE; 1135 ierr = TurnForward(ts);CHKERRQ(ierr); 1136 ierr = ReCompute(ts,tjsch,stridenum*tjsch->stride-1,stridenum*tjsch->stride);CHKERRQ(ierr); 1137 tjsch->skip_trajectory = PETSC_FALSE; 1138 } else { 1139 ierr = LoadSingle(tj,ts,stack,stridenum);CHKERRQ(ierr); 1140 ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr); 1141 tjsch->recompute = PETSC_TRUE; 1142 ierr = TurnForward(ts);CHKERRQ(ierr); 1143 ierr = ReCompute(ts,tjsch,(stridenum-1)*tjsch->stride,stridenum*tjsch->stride);CHKERRQ(ierr); 1144 } 1145 PetscFunctionReturn(0); 1146 } 1147 /* restore a checkpoint */ 1148 ierr = StackTop(stack,&e);CHKERRQ(ierr); 1149 ierr = UpdateTS(ts,stack,e);CHKERRQ(ierr); 1150 /* start with restoring a checkpoint */ 1151 tjsch->rctx->capo = stepnum; 1152 tjsch->rctx->oldcapo = tjsch->rctx->capo; 1153 shift = stepnum-localstepnum; 1154 whattodo = revolve_action(&tjsch->rctx->check,&tjsch->rctx->capo,&tjsch->rctx->fine,tjsch->rctx->snaps_in,&tjsch->rctx->info,&tjsch->rctx->where); 1155 ierr = printwhattodo(tj->monitor,whattodo,tjsch->rctx,shift);CHKERRQ(ierr); 1156 tjsch->recompute = PETSC_TRUE; 1157 ierr = TurnForward(ts);CHKERRQ(ierr); 1158 ierr = ReCompute(ts,tjsch,e->stepnum,stepnum);CHKERRQ(ierr); 1159 if (e->stepnum+1 == stepnum) { 1160 ierr = StackPop(stack,&e);CHKERRQ(ierr); 1161 ierr = ElementDestroy(stack,e);CHKERRQ(ierr); 1162 } 1163 } else { 1164 /* fill stack with info */ 1165 if (localstepnum == 0 && tjsch->total_steps-stepnum >= laststridesize) { 1166 if (tjsch->save_stack) { 1167 ierr = StackLoadAll(tj,ts,stack,stridenum);CHKERRQ(ierr); 1168 ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr); 1169 ierr = FastForwardRevolve(tjsch->rctx);CHKERRQ(ierr); 1170 } else { 1171 ierr = LoadSingle(tj,ts,stack,stridenum);CHKERRQ(ierr); 1172 ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr); 1173 ierr = ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,tjsch->total_steps,(stridenum-1)*tjsch->stride+1,1,PETSC_FALSE,&store);CHKERRQ(ierr); 1174 if (tj->monitor) { 1175 ierr = PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 1176 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); 1177 ierr = PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 1178 } 1179 ierr = ElementCreate(ts,stack,&e);CHKERRQ(ierr); 1180 ierr = ElementSet(ts,stack,&e,(stridenum-1)*tjsch->stride+1,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 1181 ierr = StackPush(stack,e);CHKERRQ(ierr); 1182 tjsch->recompute = PETSC_TRUE; 1183 ierr = TurnForward(ts);CHKERRQ(ierr); 1184 ierr = ReCompute(ts,tjsch,e->stepnum,stridenum*tjsch->stride);CHKERRQ(ierr); 1185 } 1186 PetscFunctionReturn(0); 1187 } 1188 /* restore a checkpoint */ 1189 ierr = StackTop(stack,&e);CHKERRQ(ierr); 1190 ierr = UpdateTS(ts,stack,e);CHKERRQ(ierr); 1191 /* 2 revolve actions: restore a checkpoint and then advance */ 1192 ierr = ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,tjsch->total_steps,stepnum,localstepnum,PETSC_FALSE,&store);CHKERRQ(ierr); 1193 if (tj->monitor) { 1194 ierr = PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 1195 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); 1196 ierr = PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 1197 } 1198 if (!tjsch->rctx->reverseonestep && tjsch->rctx->stepsleft > 0) tjsch->rctx->stepsleft--; 1199 if (e->stepnum < stepnum) { 1200 tjsch->recompute = PETSC_TRUE; 1201 ierr = TurnForward(ts);CHKERRQ(ierr); 1202 ierr = ReCompute(ts,tjsch,e->stepnum,stepnum);CHKERRQ(ierr); 1203 } 1204 if (e->stepnum == stepnum) { 1205 ierr = StackPop(stack,&e);CHKERRQ(ierr); 1206 ierr = ElementDestroy(stack,e);CHKERRQ(ierr); 1207 } 1208 } 1209 tjsch->rctx->reverseonestep = PETSC_FALSE; 1210 PetscFunctionReturn(0); 1211 } 1212 1213 static PetscErrorCode SetTrajTLTR(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X) 1214 { 1215 Stack *stack = &tjsch->stack; 1216 PetscInt store,localstepnum,stridenum,laststridesize; 1217 StackElement e; 1218 PetscBool done = PETSC_FALSE; 1219 PetscErrorCode ierr; 1220 1221 PetscFunctionBegin; 1222 if (!stack->solution_only && stepnum == 0) PetscFunctionReturn(0); 1223 if (stack->solution_only && stepnum == tjsch->total_steps) PetscFunctionReturn(0); 1224 1225 localstepnum = stepnum%tjsch->stride; /* index at the bottom level (inside a stride) */ 1226 stridenum = stepnum/tjsch->stride; /* index at the top level */ 1227 laststridesize = tjsch->total_steps%tjsch->stride; 1228 if (!laststridesize) laststridesize = tjsch->stride; 1229 if (stack->solution_only && localstepnum == 0 && !tjsch->rctx2->reverseonestep) { 1230 ierr = ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx2,(tjsch->total_steps+tjsch->stride-1)/tjsch->stride,stridenum,stridenum,PETSC_TRUE,&tjsch->store_stride);CHKERRQ(ierr); 1231 if (laststridesize < tjsch->stride && stepnum == tjsch->total_steps-laststridesize) { 1232 ierr = InitRevolve(laststridesize,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr); 1233 } 1234 } 1235 if (!stack->solution_only && localstepnum == 1 && !tjsch->rctx2->reverseonestep) { 1236 ierr = ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx2,(tjsch->total_steps+tjsch->stride-1)/tjsch->stride,stridenum,stridenum,PETSC_TRUE,&tjsch->store_stride);CHKERRQ(ierr); 1237 if (laststridesize < tjsch->stride && stepnum == tjsch->total_steps-laststridesize+1) { 1238 ierr = InitRevolve(laststridesize,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr); 1239 } 1240 } 1241 if (tjsch->store_stride) { 1242 ierr = TopLevelStore(tj,ts,tjsch,stepnum,localstepnum,laststridesize,&done);CHKERRQ(ierr); 1243 if (done) { 1244 ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr); 1245 PetscFunctionReturn(0); 1246 } 1247 } 1248 if (stepnum < tjsch->total_steps-laststridesize) { 1249 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 */ 1250 if (!tjsch->save_stack && !tjsch->rctx2->reverseonestep) PetscFunctionReturn(0); /* store operation does not require revolve be called at bottom level */ 1251 } 1252 /* 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 */ 1253 if (!stack->solution_only && localstepnum == 0 && stepnum != tjsch->total_steps && !tjsch->recompute) PetscFunctionReturn(0); 1254 ierr = ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,tjsch->total_steps,stepnum,localstepnum,PETSC_FALSE,&store);CHKERRQ(ierr); 1255 if (store == 1) { 1256 if (localstepnum < stack->top) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_MEMC,"Illegal modification of a non-top stack element"); 1257 ierr = ElementCreate(ts,stack,&e);CHKERRQ(ierr); 1258 ierr = ElementSet(ts,stack,&e,stepnum,time,X);CHKERRQ(ierr); 1259 ierr = StackPush(stack,e);CHKERRQ(ierr); 1260 } 1261 PetscFunctionReturn(0); 1262 } 1263 1264 static PetscErrorCode GetTrajTLTR(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum) 1265 { 1266 Stack *stack = &tjsch->stack; 1267 DiskStack *diskstack = &tjsch->diskstack; 1268 PetscInt whattodo,shift; 1269 PetscInt localstepnum,stridenum,restoredstridenum,laststridesize,store; 1270 StackElement e; 1271 PetscErrorCode ierr; 1272 1273 PetscFunctionBegin; 1274 localstepnum = stepnum%tjsch->stride; 1275 stridenum = stepnum/tjsch->stride; 1276 if (stepnum == tjsch->total_steps) { 1277 ierr = TurnBackward(ts);CHKERRQ(ierr); 1278 tjsch->rctx->reverseonestep = PETSC_FALSE; 1279 PetscFunctionReturn(0); 1280 } 1281 laststridesize = tjsch->total_steps%tjsch->stride; 1282 if (!laststridesize) laststridesize = tjsch->stride; 1283 /* 1284 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: 1285 Case 1 (save_stack) 1286 Restore a disk checkpoint; update TS with the last element in the restored data; recompute to the current point. 1287 Case 2 (!save_stack) 1288 Restore a disk checkpoint; update TS with the restored point; recompute to the current point. 1289 */ 1290 if (localstepnum == 0 && stepnum <= tjsch->total_steps-laststridesize) { 1291 /* restore the top element in the stack for disk checkpoints */ 1292 restoredstridenum = diskstack->container[diskstack->top]; 1293 tjsch->rctx2->reverseonestep = PETSC_FALSE; 1294 /* top-level revolve must be applied before current step, just like the solution_only mode for single-level revolve */ 1295 if (!tjsch->save_stack && stack->solution_only) { /* start with restoring a checkpoint */ 1296 tjsch->rctx2->capo = stridenum; 1297 tjsch->rctx2->oldcapo = tjsch->rctx2->capo; 1298 shift = 0; 1299 whattodo = revolve2_action(&tjsch->rctx2->check,&tjsch->rctx2->capo,&tjsch->rctx2->fine,tjsch->rctx2->snaps_in,&tjsch->rctx2->info,&tjsch->rctx2->where); 1300 ierr = printwhattodo2(tj->monitor,whattodo,tjsch->rctx2,shift);CHKERRQ(ierr); 1301 } else { /* 2 revolve actions: restore a checkpoint and then advance */ 1302 ierr = ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx2,(tjsch->total_steps+tjsch->stride-1)/tjsch->stride,stridenum,stridenum,PETSC_TRUE,&tjsch->store_stride);CHKERRQ(ierr); 1303 if (tj->monitor) { 1304 ierr = PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 1305 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); 1306 ierr = PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 1307 } 1308 if (!tjsch->rctx2->reverseonestep && tjsch->rctx2->stepsleft > 0) tjsch->rctx2->stepsleft--; 1309 } 1310 /* fill stack */ 1311 if (stack->solution_only) { 1312 if (tjsch->save_stack) { 1313 if (restoredstridenum < stridenum) { 1314 ierr = StackLoadLast(tj,ts,stack,restoredstridenum);CHKERRQ(ierr); 1315 } else { 1316 ierr = StackLoadAll(tj,ts,stack,restoredstridenum);CHKERRQ(ierr); 1317 } 1318 /* recompute one step ahead */ 1319 tjsch->recompute = PETSC_TRUE; 1320 tjsch->skip_trajectory = PETSC_TRUE; 1321 ierr = TurnForward(ts);CHKERRQ(ierr); 1322 ierr = ReCompute(ts,tjsch,stridenum*tjsch->stride-1,stridenum*tjsch->stride);CHKERRQ(ierr); 1323 tjsch->skip_trajectory = PETSC_FALSE; 1324 if (restoredstridenum < stridenum) { 1325 ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr); 1326 tjsch->recompute = PETSC_TRUE; 1327 ierr = TurnForward(ts);CHKERRQ(ierr); 1328 ierr = ReCompute(ts,tjsch,restoredstridenum*tjsch->stride,stepnum);CHKERRQ(ierr); 1329 } else { /* stack ready, fast forward revolve status */ 1330 ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr); 1331 ierr = FastForwardRevolve(tjsch->rctx);CHKERRQ(ierr); 1332 } 1333 } else { 1334 ierr = LoadSingle(tj,ts,stack,restoredstridenum);CHKERRQ(ierr); 1335 ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr); 1336 tjsch->recompute = PETSC_TRUE; 1337 ierr = TurnForward(ts);CHKERRQ(ierr); 1338 ierr = ReCompute(ts,tjsch,(restoredstridenum-1)*tjsch->stride,stepnum);CHKERRQ(ierr); 1339 } 1340 } else { 1341 if (tjsch->save_stack) { 1342 if (restoredstridenum < stridenum) { 1343 ierr = StackLoadLast(tj,ts,stack,restoredstridenum);CHKERRQ(ierr); 1344 /* reset revolve */ 1345 ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr); 1346 tjsch->recompute = PETSC_TRUE; 1347 ierr = TurnForward(ts);CHKERRQ(ierr); 1348 ierr = ReCompute(ts,tjsch,restoredstridenum*tjsch->stride,stepnum);CHKERRQ(ierr); 1349 } else { /* stack ready, fast forward revolve status */ 1350 ierr = StackLoadAll(tj,ts,stack,restoredstridenum);CHKERRQ(ierr); 1351 ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr); 1352 ierr = FastForwardRevolve(tjsch->rctx);CHKERRQ(ierr); 1353 } 1354 } else { 1355 ierr = LoadSingle(tj,ts,stack,restoredstridenum);CHKERRQ(ierr); 1356 ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr); 1357 /* push first element to stack */ 1358 if (tjsch->store_stride || tjsch->rctx2->reverseonestep) { 1359 shift = (restoredstridenum-1)*tjsch->stride-localstepnum; 1360 ierr = ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,tjsch->total_steps,(restoredstridenum-1)*tjsch->stride+1,1,PETSC_FALSE,&store);CHKERRQ(ierr); 1361 if (tj->monitor) { 1362 ierr = PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 1363 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); 1364 ierr = PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 1365 } 1366 ierr = ElementCreate(ts,stack,&e);CHKERRQ(ierr); 1367 ierr = ElementSet(ts,stack,&e,(restoredstridenum-1)*tjsch->stride+1,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 1368 ierr = StackPush(stack,e);CHKERRQ(ierr); 1369 } 1370 tjsch->recompute = PETSC_TRUE; 1371 ierr = TurnForward(ts);CHKERRQ(ierr); 1372 ierr = ReCompute(ts,tjsch,(restoredstridenum-1)*tjsch->stride+1,stepnum);CHKERRQ(ierr); 1373 } 1374 } 1375 if (restoredstridenum == stridenum) diskstack->top--; 1376 tjsch->rctx->reverseonestep = PETSC_FALSE; 1377 PetscFunctionReturn(0); 1378 } 1379 1380 if (stack->solution_only) { 1381 /* restore a checkpoint */ 1382 ierr = StackTop(stack,&e);CHKERRQ(ierr); 1383 ierr = UpdateTS(ts,stack,e);CHKERRQ(ierr); 1384 /* start with restoring a checkpoint */ 1385 tjsch->rctx->capo = stepnum; 1386 tjsch->rctx->oldcapo = tjsch->rctx->capo; 1387 shift = stepnum-localstepnum; 1388 whattodo = revolve_action(&tjsch->rctx->check,&tjsch->rctx->capo,&tjsch->rctx->fine,tjsch->rctx->snaps_in,&tjsch->rctx->info,&tjsch->rctx->where); 1389 ierr = printwhattodo(tj->monitor,whattodo,tjsch->rctx,shift);CHKERRQ(ierr); 1390 tjsch->recompute = PETSC_TRUE; 1391 ierr = TurnForward(ts);CHKERRQ(ierr); 1392 ierr = ReCompute(ts,tjsch,e->stepnum,stepnum);CHKERRQ(ierr); 1393 if (e->stepnum+1 == stepnum) { 1394 ierr = StackPop(stack,&e);CHKERRQ(ierr); 1395 ierr = ElementDestroy(stack,e);CHKERRQ(ierr); 1396 } 1397 } else { 1398 /* restore a checkpoint */ 1399 ierr = StackTop(stack,&e);CHKERRQ(ierr); 1400 ierr = UpdateTS(ts,stack,e);CHKERRQ(ierr); 1401 /* 2 revolve actions: restore a checkpoint and then advance */ 1402 ierr = ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,tjsch->total_steps,stepnum,localstepnum,PETSC_FALSE,&store);CHKERRQ(ierr); 1403 if (tj->monitor) { 1404 ierr = PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 1405 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); 1406 ierr = PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 1407 } 1408 if (!tjsch->rctx->reverseonestep && tjsch->rctx->stepsleft > 0) tjsch->rctx->stepsleft--; 1409 if (e->stepnum < stepnum) { 1410 tjsch->recompute = PETSC_TRUE; 1411 ierr = TurnForward(ts);CHKERRQ(ierr); 1412 ierr = ReCompute(ts,tjsch,e->stepnum,stepnum);CHKERRQ(ierr); 1413 } 1414 if (e->stepnum == stepnum) { 1415 ierr = StackPop(stack,&e);CHKERRQ(ierr); 1416 ierr = ElementDestroy(stack,e);CHKERRQ(ierr); 1417 } 1418 } 1419 tjsch->rctx->reverseonestep = PETSC_FALSE; 1420 PetscFunctionReturn(0); 1421 } 1422 1423 static PetscErrorCode SetTrajRMS(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X) 1424 { 1425 Stack *stack = &tjsch->stack; 1426 PetscInt store; 1427 StackElement e; 1428 PetscErrorCode ierr; 1429 1430 PetscFunctionBegin; 1431 if (!stack->solution_only && stepnum == 0) PetscFunctionReturn(0); 1432 if (stack->solution_only && stepnum == tjsch->total_steps) PetscFunctionReturn(0); 1433 ierr = ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,tjsch->total_steps,stepnum,stepnum,PETSC_FALSE,&store);CHKERRQ(ierr); 1434 if (store == 1){ 1435 if (stepnum < stack->top) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_MEMC,"Illegal modification of a non-top stack element"); 1436 ierr = ElementCreate(ts,stack,&e);CHKERRQ(ierr); 1437 ierr = ElementSet(ts,stack,&e,stepnum,time,X);CHKERRQ(ierr); 1438 ierr = StackPush(stack,e);CHKERRQ(ierr); 1439 } else if (store == 2) { 1440 ierr = DumpSingle(tj,ts,stack,tjsch->rctx->check+1);CHKERRQ(ierr); 1441 } 1442 PetscFunctionReturn(0); 1443 } 1444 1445 static PetscErrorCode GetTrajRMS(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum) 1446 { 1447 Stack *stack = &tjsch->stack; 1448 PetscInt whattodo,shift; 1449 PetscInt restart; 1450 PetscBool ondisk; 1451 StackElement e; 1452 PetscErrorCode ierr; 1453 1454 PetscFunctionBegin; 1455 if (stepnum == 0 || stepnum == tjsch->total_steps) { 1456 ierr = TurnBackward(ts);CHKERRQ(ierr); 1457 tjsch->rctx->reverseonestep = PETSC_FALSE; 1458 PetscFunctionReturn(0); 1459 } 1460 tjsch->rctx->capo = stepnum; 1461 tjsch->rctx->oldcapo = tjsch->rctx->capo; 1462 shift = 0; 1463 whattodo = revolve_action(&tjsch->rctx->check,&tjsch->rctx->capo,&tjsch->rctx->fine,tjsch->rctx->snaps_in,&tjsch->rctx->info,&tjsch->rctx->where); /* whattodo=restore */ 1464 ierr = printwhattodo(tj->monitor,whattodo,tjsch->rctx,shift);CHKERRQ(ierr); 1465 /* restore a checkpoint */ 1466 restart = tjsch->rctx->capo; 1467 if (!tjsch->rctx->where) { 1468 ondisk = PETSC_TRUE; 1469 ierr = LoadSingle(tj,ts,stack,tjsch->rctx->check+1);CHKERRQ(ierr); 1470 ierr = TurnBackward(ts);CHKERRQ(ierr); 1471 } else { 1472 ondisk = PETSC_FALSE; 1473 ierr = StackTop(stack,&e);CHKERRQ(ierr); 1474 ierr = UpdateTS(ts,stack,e);CHKERRQ(ierr); 1475 } 1476 if (!stack->solution_only) { /* whattodo must be 5 or 8 */ 1477 /* ask Revolve what to do next */ 1478 tjsch->rctx->oldcapo = tjsch->rctx->capo; 1479 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*/ 1480 ierr = printwhattodo(tj->monitor,whattodo,tjsch->rctx,shift);CHKERRQ(ierr); 1481 if (whattodo == 3 || whattodo == 4) tjsch->rctx->reverseonestep = PETSC_TRUE; 1482 if (whattodo == 1) tjsch->rctx->stepsleft = tjsch->rctx->capo-tjsch->rctx->oldcapo; 1483 if (tj->monitor) { 1484 ierr = PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 1485 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); 1486 ierr = PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 1487 } 1488 if (!tjsch->rctx->reverseonestep && tjsch->rctx->stepsleft > 0) tjsch->rctx->stepsleft--; 1489 restart++; /* skip one step */ 1490 } 1491 if (stack->solution_only || (!stack->solution_only && restart < stepnum)) { 1492 tjsch->recompute = PETSC_TRUE; 1493 ierr = TurnForward(ts);CHKERRQ(ierr); 1494 ierr = ReCompute(ts,tjsch,restart,stepnum);CHKERRQ(ierr); 1495 } 1496 if (!ondisk && ( (stack->solution_only && e->stepnum+1 == stepnum) || (!stack->solution_only && e->stepnum == stepnum) )) { 1497 ierr = StackPop(stack,&e);CHKERRQ(ierr); 1498 ierr = ElementDestroy(stack,e);CHKERRQ(ierr); 1499 } 1500 tjsch->rctx->reverseonestep = PETSC_FALSE; 1501 PetscFunctionReturn(0); 1502 } 1503 #endif 1504 1505 static PetscErrorCode TSTrajectorySet_Memory(TSTrajectory tj,TS ts,PetscInt stepnum,PetscReal time,Vec X) 1506 { 1507 TJScheduler *tjsch = (TJScheduler*)tj->data; 1508 PetscErrorCode ierr; 1509 1510 PetscFunctionBegin; 1511 if (!tjsch->recompute) { /* use global stepnum in the forward sweep */ 1512 ierr = TSGetTotalSteps(ts,&stepnum);CHKERRQ(ierr); 1513 } 1514 /* for consistency */ 1515 if (!tjsch->recompute && stepnum == 0) ts->ptime_prev = ts->ptime-ts->time_step; 1516 switch (tjsch->stype) { 1517 case NONE: 1518 ierr = SetTrajN(ts,tjsch,stepnum,time,X);CHKERRQ(ierr); 1519 break; 1520 case TWO_LEVEL_NOREVOLVE: 1521 ierr = SetTrajTLNR(tj,ts,tjsch,stepnum,time,X);CHKERRQ(ierr); 1522 break; 1523 #if defined(PETSC_HAVE_REVOLVE) 1524 case TWO_LEVEL_REVOLVE: 1525 ierr = SetTrajTLR(tj,ts,tjsch,stepnum,time,X);CHKERRQ(ierr); 1526 break; 1527 case TWO_LEVEL_TWO_REVOLVE: 1528 ierr = SetTrajTLTR(tj,ts,tjsch,stepnum,time,X);CHKERRQ(ierr); 1529 break; 1530 case REVOLVE_OFFLINE: 1531 ierr = SetTrajROF(tj,ts,tjsch,stepnum,time,X);CHKERRQ(ierr); 1532 break; 1533 case REVOLVE_ONLINE: 1534 ierr = SetTrajRON(tj,ts,tjsch,stepnum,time,X);CHKERRQ(ierr); 1535 break; 1536 case REVOLVE_MULTISTAGE: 1537 ierr = SetTrajRMS(tj,ts,tjsch,stepnum,time,X);CHKERRQ(ierr); 1538 break; 1539 #endif 1540 default: 1541 break; 1542 } 1543 PetscFunctionReturn(0); 1544 } 1545 1546 static PetscErrorCode TSTrajectoryGet_Memory(TSTrajectory tj,TS ts,PetscInt stepnum,PetscReal *t) 1547 { 1548 TJScheduler *tjsch = (TJScheduler*)tj->data; 1549 PetscErrorCode ierr; 1550 1551 PetscFunctionBegin; 1552 ierr = TSGetTotalSteps(ts,&stepnum);CHKERRQ(ierr); 1553 if (stepnum == 0) PetscFunctionReturn(0); 1554 switch (tjsch->stype) { 1555 case NONE: 1556 ierr = GetTrajN(ts,tjsch,stepnum);CHKERRQ(ierr); 1557 break; 1558 case TWO_LEVEL_NOREVOLVE: 1559 ierr = GetTrajTLNR(tj,ts,tjsch,stepnum);CHKERRQ(ierr); 1560 break; 1561 #if defined(PETSC_HAVE_REVOLVE) 1562 case TWO_LEVEL_REVOLVE: 1563 ierr = GetTrajTLR(tj,ts,tjsch,stepnum);CHKERRQ(ierr); 1564 break; 1565 case TWO_LEVEL_TWO_REVOLVE: 1566 ierr = GetTrajTLTR(tj,ts,tjsch,stepnum);CHKERRQ(ierr); 1567 break; 1568 case REVOLVE_OFFLINE: 1569 ierr = GetTrajROF(tj,ts,tjsch,stepnum);CHKERRQ(ierr); 1570 break; 1571 case REVOLVE_ONLINE: 1572 ierr = GetTrajRON(tj,ts,tjsch,stepnum);CHKERRQ(ierr); 1573 break; 1574 case REVOLVE_MULTISTAGE: 1575 ierr = GetTrajRMS(tj,ts,tjsch,stepnum);CHKERRQ(ierr); 1576 break; 1577 #endif 1578 default: 1579 break; 1580 } 1581 PetscFunctionReturn(0); 1582 } 1583 1584 PETSC_UNUSED static PetscErrorCode TSTrajectorySetStride_Memory(TSTrajectory tj,PetscInt stride) 1585 { 1586 TJScheduler *tjsch = (TJScheduler*)tj->data; 1587 1588 PetscFunctionBegin; 1589 tjsch->stride = stride; 1590 PetscFunctionReturn(0); 1591 } 1592 1593 PETSC_UNUSED static PetscErrorCode TSTrajectorySetMaxCpsRAM_Memory(TSTrajectory tj,PetscInt max_cps_ram) 1594 { 1595 TJScheduler *tjsch = (TJScheduler*)tj->data; 1596 1597 PetscFunctionBegin; 1598 tjsch->max_cps_ram = max_cps_ram; 1599 PetscFunctionReturn(0); 1600 } 1601 1602 PETSC_UNUSED static PetscErrorCode TSTrajectorySetMaxCpsDisk_Memory(TSTrajectory tj,PetscInt max_cps_disk) 1603 { 1604 TJScheduler *tjsch = (TJScheduler*)tj->data; 1605 1606 PetscFunctionBegin; 1607 tjsch->max_cps_disk = max_cps_disk; 1608 PetscFunctionReturn(0); 1609 } 1610 1611 #if defined(PETSC_HAVE_REVOLVE) 1612 PETSC_UNUSED static PetscErrorCode TSTrajectorySetRevolveOnline(TSTrajectory tj,PetscBool use_online) 1613 { 1614 TJScheduler *tjsch = (TJScheduler*)tj->data; 1615 1616 PetscFunctionBegin; 1617 tjsch->use_online = use_online; 1618 PetscFunctionReturn(0); 1619 } 1620 #endif 1621 1622 PETSC_UNUSED static PetscErrorCode TSTrajectorySetSaveStack(TSTrajectory tj,PetscBool save_stack) 1623 { 1624 TJScheduler *tjsch = (TJScheduler*)tj->data; 1625 1626 PetscFunctionBegin; 1627 tjsch->save_stack = save_stack; 1628 PetscFunctionReturn(0); 1629 } 1630 1631 PETSC_UNUSED static PetscErrorCode TSTrajectorySetSolutionOnly(TSTrajectory tj,PetscBool solution_only) 1632 { 1633 TJScheduler *tjsch = (TJScheduler*)tj->data; 1634 Stack *stack = &tjsch->stack; 1635 1636 PetscFunctionBegin; 1637 stack->solution_only = solution_only; 1638 PetscFunctionReturn(0); 1639 } 1640 1641 PETSC_UNUSED static PetscErrorCode TSTrajectorySetUseDRAM(TSTrajectory tj,PetscBool use_dram) 1642 { 1643 TJScheduler *tjsch = (TJScheduler*)tj->data; 1644 1645 PetscFunctionBegin; 1646 tjsch->stack.use_dram = use_dram; 1647 PetscFunctionReturn(0); 1648 } 1649 1650 static PetscErrorCode TSTrajectorySetFromOptions_Memory(PetscOptionItems *PetscOptionsObject,TSTrajectory tj) 1651 { 1652 TJScheduler *tjsch = (TJScheduler*)tj->data; 1653 PetscErrorCode ierr; 1654 1655 PetscFunctionBegin; 1656 ierr = PetscOptionsHead(PetscOptionsObject,"Memory based TS trajectory options");CHKERRQ(ierr); 1657 { 1658 ierr = PetscOptionsInt("-ts_trajectory_max_cps_ram","Maximum number of checkpoints in RAM","TSTrajectorySetMaxCpsRAM_Memory",tjsch->max_cps_ram,&tjsch->max_cps_ram,NULL);CHKERRQ(ierr); 1659 ierr = PetscOptionsInt("-ts_trajectory_max_cps_disk","Maximum number of checkpoints on disk","TSTrajectorySetMaxCpsDisk_Memory",tjsch->max_cps_disk,&tjsch->max_cps_disk,NULL);CHKERRQ(ierr); 1660 ierr = PetscOptionsInt("-ts_trajectory_stride","Stride to save checkpoints to file","TSTrajectorySetStride_Memory",tjsch->stride,&tjsch->stride,NULL);CHKERRQ(ierr); 1661 #if defined(PETSC_HAVE_REVOLVE) 1662 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); 1663 #endif 1664 ierr = PetscOptionsBool("-ts_trajectory_save_stack","Save all stack to disk","TSTrajectorySetSaveStack",tjsch->save_stack,&tjsch->save_stack,NULL);CHKERRQ(ierr); 1665 ierr = PetscOptionsBool("-ts_trajectory_solution_only","Checkpoint solution only","TSTrajectorySetSolutionOnly",tjsch->stack.solution_only,&tjsch->stack.solution_only,NULL);CHKERRQ(ierr); 1666 ierr = PetscOptionsBool("-ts_trajectory_use_dram","Use DRAM for checkpointing","TSTrajectorySetUseDRAM",tjsch->stack.use_dram,&tjsch->stack.use_dram,NULL);CHKERRQ(ierr); 1667 } 1668 ierr = PetscOptionsTail();CHKERRQ(ierr); 1669 PetscFunctionReturn(0); 1670 } 1671 1672 static PetscErrorCode TSTrajectorySetUp_Memory(TSTrajectory tj,TS ts) 1673 { 1674 TJScheduler *tjsch = (TJScheduler*)tj->data; 1675 Stack *stack = &tjsch->stack; 1676 #if defined(PETSC_HAVE_REVOLVE) 1677 RevolveCTX *rctx,*rctx2; 1678 DiskStack *diskstack = &tjsch->diskstack; 1679 PetscInt diskblocks; 1680 #endif 1681 PetscInt numY; 1682 PetscBool flg; 1683 PetscErrorCode ierr; 1684 1685 PetscFunctionBegin; 1686 PetscStrcmp(((PetscObject)ts->adapt)->type_name,TSADAPTNONE,&flg); 1687 if (flg) tjsch->total_steps = PetscMin(ts->max_steps,(PetscInt)(ceil(ts->max_time/ts->time_step))); /* fixed time step */ 1688 if (tjsch->max_cps_ram > 0) stack->stacksize = tjsch->max_cps_ram; 1689 1690 if (tjsch->stride > 1) { /* two level mode */ 1691 if (tjsch->save_stack && tjsch->max_cps_disk > 1 && tjsch->max_cps_disk <= tjsch->max_cps_ram) SETERRQ(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."); 1692 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 */ 1693 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 */ 1694 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 */ 1695 } else { /* single level mode */ 1696 if (flg) { /* fixed time step */ 1697 if (tjsch->max_cps_ram >= tjsch->total_steps-1 || tjsch->max_cps_ram < 1) tjsch->stype = NONE; /* checkpoint all */ 1698 else tjsch->stype = (tjsch->max_cps_disk>1) ? REVOLVE_MULTISTAGE : REVOLVE_OFFLINE; 1699 } else tjsch->stype = NONE; /* checkpoint all for adaptive time step */ 1700 #if defined(PETSC_HAVE_REVOLVE) 1701 if (tjsch->use_online) tjsch->stype = REVOLVE_ONLINE; /* trick into online (for testing purpose only) */ 1702 #endif 1703 } 1704 1705 if (tjsch->stype > TWO_LEVEL_NOREVOLVE) { 1706 #ifndef PETSC_HAVE_REVOLVE 1707 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."); 1708 #else 1709 switch (tjsch->stype) { 1710 case TWO_LEVEL_REVOLVE: 1711 revolve_create_offline(tjsch->stride,tjsch->max_cps_ram); 1712 break; 1713 case TWO_LEVEL_TWO_REVOLVE: 1714 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. */ 1715 diskstack->stacksize = diskblocks; 1716 revolve_create_offline(tjsch->stride,tjsch->max_cps_ram); 1717 revolve2_create_offline((tjsch->total_steps+tjsch->stride-1)/tjsch->stride,diskblocks); 1718 ierr = PetscCalloc1(1,&rctx2);CHKERRQ(ierr); 1719 rctx2->snaps_in = diskblocks; 1720 rctx2->reverseonestep = PETSC_FALSE; 1721 rctx2->check = 0; 1722 rctx2->oldcapo = 0; 1723 rctx2->capo = 0; 1724 rctx2->info = 2; 1725 rctx2->fine = (tjsch->total_steps+tjsch->stride-1)/tjsch->stride; 1726 tjsch->rctx2 = rctx2; 1727 diskstack->top = -1; 1728 ierr = PetscMalloc1(diskstack->stacksize*sizeof(PetscInt),&diskstack->container);CHKERRQ(ierr); 1729 break; 1730 case REVOLVE_OFFLINE: 1731 revolve_create_offline(tjsch->total_steps,tjsch->max_cps_ram); 1732 break; 1733 case REVOLVE_ONLINE: 1734 stack->stacksize = tjsch->max_cps_ram; 1735 revolve_create_online(tjsch->max_cps_ram); 1736 break; 1737 case REVOLVE_MULTISTAGE: 1738 revolve_create_multistage(tjsch->total_steps,tjsch->max_cps_ram+tjsch->max_cps_disk,tjsch->max_cps_ram); 1739 break; 1740 default: 1741 break; 1742 } 1743 ierr = PetscCalloc1(1,&rctx);CHKERRQ(ierr); 1744 rctx->snaps_in = tjsch->max_cps_ram; /* for theta methods snaps_in=2*max_cps_ram */ 1745 rctx->reverseonestep = PETSC_FALSE; 1746 rctx->check = 0; 1747 rctx->oldcapo = 0; 1748 rctx->capo = 0; 1749 rctx->info = 2; 1750 rctx->fine = (tjsch->stride > 1) ? tjsch->stride : tjsch->total_steps; 1751 tjsch->rctx = rctx; 1752 if (tjsch->stype == REVOLVE_ONLINE) rctx->fine = -1; 1753 #endif 1754 } else { 1755 if (tjsch->stype == TWO_LEVEL_NOREVOLVE) stack->stacksize = tjsch->stride-1; /* need tjsch->stride-1 at most */ 1756 if (tjsch->stype == NONE) { 1757 if (flg) stack->stacksize = stack->solution_only ? tjsch->total_steps : tjsch->total_steps-1; /* fix time step */ 1758 else { /* adaptive time step */ 1759 if(tjsch->max_cps_ram == -1) stack->stacksize = ts->max_steps; /* if max_cps_ram is not specified, use maximal allowed number of steps for stack size */ 1760 tjsch->total_steps = stack->solution_only ? stack->stacksize:stack->stacksize+1; /* will be updated as time integration advances */ 1761 } 1762 } 1763 } 1764 1765 tjsch->recompute = PETSC_FALSE; 1766 ierr = TSGetStages(ts,&numY,PETSC_IGNORE);CHKERRQ(ierr); 1767 ierr = StackCreate(stack,stack->stacksize,numY);CHKERRQ(ierr); 1768 PetscFunctionReturn(0); 1769 } 1770 1771 static PetscErrorCode TSTrajectoryDestroy_Memory(TSTrajectory tj) 1772 { 1773 TJScheduler *tjsch = (TJScheduler*)tj->data; 1774 PetscErrorCode ierr; 1775 1776 PetscFunctionBegin; 1777 if (tjsch->stype > TWO_LEVEL_NOREVOLVE) { 1778 #if defined(PETSC_HAVE_REVOLVE) 1779 revolve_reset(); 1780 if (tjsch->stype == TWO_LEVEL_TWO_REVOLVE) { 1781 revolve2_reset(); 1782 ierr = PetscFree(tjsch->diskstack.container);CHKERRQ(ierr); 1783 } 1784 #endif 1785 } 1786 ierr = StackDestroy(&tjsch->stack);CHKERRQ(ierr); 1787 #if defined(PETSC_HAVE_REVOLVE) 1788 if (tjsch->stype > TWO_LEVEL_NOREVOLVE) { 1789 ierr = PetscFree(tjsch->rctx);CHKERRQ(ierr); 1790 ierr = PetscFree(tjsch->rctx2);CHKERRQ(ierr); 1791 } 1792 #endif 1793 ierr = PetscFree(tjsch);CHKERRQ(ierr); 1794 PetscFunctionReturn(0); 1795 } 1796 1797 /*MC 1798 TSTRAJECTORYMEMORY - Stores each solution of the ODE/ADE in memory 1799 1800 Level: intermediate 1801 1802 .seealso: TSTrajectoryCreate(), TS, TSTrajectorySetType() 1803 1804 M*/ 1805 PETSC_EXTERN PetscErrorCode TSTrajectoryCreate_Memory(TSTrajectory tj,TS ts) 1806 { 1807 TJScheduler *tjsch; 1808 PetscErrorCode ierr; 1809 1810 PetscFunctionBegin; 1811 tj->ops->set = TSTrajectorySet_Memory; 1812 tj->ops->get = TSTrajectoryGet_Memory; 1813 tj->ops->setup = TSTrajectorySetUp_Memory; 1814 tj->ops->destroy = TSTrajectoryDestroy_Memory; 1815 tj->ops->setfromoptions = TSTrajectorySetFromOptions_Memory; 1816 1817 ierr = PetscCalloc1(1,&tjsch);CHKERRQ(ierr); 1818 tjsch->stype = NONE; 1819 tjsch->max_cps_ram = -1; /* -1 indicates that it is not set */ 1820 tjsch->max_cps_disk = -1; /* -1 indicates that it is not set */ 1821 tjsch->stride = 0; /* if not zero, two-level checkpointing will be used */ 1822 #if defined(PETSC_HAVE_REVOLVE) 1823 tjsch->use_online = PETSC_FALSE; 1824 #endif 1825 tjsch->save_stack = PETSC_TRUE; 1826 1827 tjsch->stack.solution_only = PETSC_TRUE; 1828 1829 tj->data = tjsch; 1830 PetscFunctionReturn(0); 1831 } 1832