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->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->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->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 = TSGetStepNumber(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->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; 518 PetscErrorCode ierr; 519 520 PetscFunctionBegin; 521 ts->steps = stepnumbegin; /* global step number */ 522 for (i=stepnumbegin;i<stepnumend;i++) { /* assume fixed step size */ 523 if (stack->solution_only && !tjsch->skip_trajectory) { /* revolve online need this */ 524 ierr = TSTrajectorySet(ts->trajectory,ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 525 } 526 ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 527 ierr = TSStep(ts);CHKERRQ(ierr); 528 if (!stack->solution_only && !tjsch->skip_trajectory) { 529 ierr = TSTrajectorySet(ts->trajectory,ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 530 } 531 ierr = TSEventHandler(ts);CHKERRQ(ierr); 532 if (!ts->steprollback) { 533 ierr = TSPostStep(ts);CHKERRQ(ierr); 534 } 535 } 536 ierr = TurnBackward(ts);CHKERRQ(ierr); 537 ts->trajectory->recomps += stepnumend-stepnumbegin; /* recomputation counter */ 538 ts->steps = stepnumend; 539 PetscFunctionReturn(0); 540 } 541 542 static PetscErrorCode TopLevelStore(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscInt localstepnum,PetscInt laststridesize,PetscBool *done) 543 { 544 Stack *stack = &tjsch->stack; 545 DiskStack *diskstack = &tjsch->diskstack; 546 PetscInt stridenum; 547 PetscErrorCode ierr; 548 549 PetscFunctionBegin; 550 *done = PETSC_FALSE; 551 stridenum = stepnum/tjsch->stride; 552 /* make sure saved checkpoint id starts from 1 553 skip last stride when using stridenum+1 554 skip first stride when using stridenum */ 555 if (stack->solution_only) { 556 if (tjsch->save_stack) { 557 if (localstepnum == tjsch->stride-1 && stepnum < tjsch->total_steps-laststridesize) { /* current step will be saved without going through stack */ 558 ierr = StackDumpAll(tj,ts,stack,stridenum+1);CHKERRQ(ierr); 559 if (tjsch->stype == TWO_LEVEL_TWO_REVOLVE) diskstack->container[++diskstack->top] = stridenum+1; 560 *done = PETSC_TRUE; 561 } 562 } else { 563 if (localstepnum == 0 && stepnum < tjsch->total_steps-laststridesize) { 564 ierr = DumpSingle(tj,ts,stack,stridenum+1);CHKERRQ(ierr); 565 if (tjsch->stype == TWO_LEVEL_TWO_REVOLVE) diskstack->container[++diskstack->top] = stridenum+1; 566 *done = PETSC_TRUE; 567 } 568 } 569 } else { 570 if (tjsch->save_stack) { 571 if (localstepnum == 0 && stepnum < tjsch->total_steps && stepnum != 0) { /* skip the first stride */ 572 ierr = StackDumpAll(tj,ts,stack,stridenum);CHKERRQ(ierr); 573 if (tjsch->stype == TWO_LEVEL_TWO_REVOLVE) diskstack->container[++diskstack->top] = stridenum; 574 *done = PETSC_TRUE; 575 } 576 } else { 577 if (localstepnum == 1 && stepnum < tjsch->total_steps-laststridesize) { 578 ierr = DumpSingle(tj,ts,stack,stridenum+1);CHKERRQ(ierr); 579 if (tjsch->stype == TWO_LEVEL_TWO_REVOLVE) diskstack->container[++diskstack->top] = stridenum+1; 580 *done = PETSC_TRUE; 581 } 582 } 583 } 584 PetscFunctionReturn(0); 585 } 586 587 static PetscErrorCode SetTrajN(TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X) 588 { 589 Stack *stack = &tjsch->stack; 590 StackElement e; 591 PetscErrorCode ierr; 592 593 PetscFunctionBegin; 594 /* skip the last step */ 595 if (ts->reason) { /* only affect the forward run */ 596 /* update total_steps in the end of forward run */ 597 if (stepnum != tjsch->total_steps) tjsch->total_steps = stepnum; 598 if (stack->solution_only) { 599 /* get rid of the solution at second last step */ 600 ierr = StackPop(stack,&e);CHKERRQ(ierr); 601 ierr = ElementDestroy(stack,e);CHKERRQ(ierr); 602 } 603 PetscFunctionReturn(0); 604 } 605 /* do not save trajectory at the recompute stage for solution_only mode */ 606 if (stack->solution_only && tjsch->recompute) PetscFunctionReturn(0); 607 /* skip the first step for no_solution_only mode */ 608 if (!stack->solution_only && stepnum == 0) PetscFunctionReturn(0); 609 610 /* resize the stack */ 611 if (stack->top+1 == stack->stacksize) { 612 ierr = StackResize(stack,2*stack->stacksize);CHKERRQ(ierr); 613 } 614 /* update timenext for the previous step; necessary for step adaptivity */ 615 if (stack->top > -1) { 616 ierr = StackTop(stack,&e);CHKERRQ(ierr); 617 e->timenext = ts->ptime; 618 } 619 if (stepnum < stack->top) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_MEMC,"Illegal modification of a non-top stack element"); 620 ierr = ElementCreate(ts,stack,&e);CHKERRQ(ierr); 621 ierr = ElementSet(ts,stack,&e,stepnum,time,X);CHKERRQ(ierr); 622 ierr = StackPush(stack,e);CHKERRQ(ierr); 623 PetscFunctionReturn(0); 624 } 625 626 static PetscErrorCode GetTrajN(TS ts,TJScheduler *tjsch,PetscInt stepnum) 627 { 628 Stack *stack = &tjsch->stack; 629 StackElement e; 630 PetscErrorCode ierr; 631 632 PetscFunctionBegin; 633 if (stepnum == tjsch->total_steps) { 634 ierr = TurnBackward(ts);CHKERRQ(ierr); 635 PetscFunctionReturn(0); 636 } 637 /* restore a checkpoint */ 638 ierr = StackTop(stack,&e);CHKERRQ(ierr); 639 ierr = UpdateTS(ts,stack,e);CHKERRQ(ierr); 640 if (stack->solution_only) {/* recompute one step */ 641 tjsch->recompute = PETSC_TRUE; 642 ierr = TurnForwardWithStepsize(ts,e->timenext-e->time);CHKERRQ(ierr); 643 ierr = ReCompute(ts,tjsch,e->stepnum,stepnum);CHKERRQ(ierr); 644 } 645 ierr = StackPop(stack,&e);CHKERRQ(ierr); 646 ierr = ElementDestroy(stack,e);CHKERRQ(ierr); 647 PetscFunctionReturn(0); 648 } 649 650 static PetscErrorCode SetTrajTLNR(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X) 651 { 652 Stack *stack = &tjsch->stack; 653 PetscInt localstepnum,laststridesize; 654 StackElement e; 655 PetscBool done; 656 PetscErrorCode ierr; 657 658 PetscFunctionBegin; 659 if (!stack->solution_only && stepnum == 0) PetscFunctionReturn(0); 660 if (stack->solution_only && stepnum == tjsch->total_steps) PetscFunctionReturn(0); 661 if (tjsch->save_stack && tjsch->recompute) PetscFunctionReturn(0); 662 663 localstepnum = stepnum%tjsch->stride; 664 /* (stride size-1) checkpoints are saved in each stride; an extra point is added by StackDumpAll() */ 665 laststridesize = tjsch->total_steps%tjsch->stride; 666 if (!laststridesize) laststridesize = tjsch->stride; 667 668 if (!tjsch->recompute) { 669 ierr = TopLevelStore(tj,ts,tjsch,stepnum,localstepnum,laststridesize,&done);CHKERRQ(ierr); 670 if (!tjsch->save_stack && stepnum < tjsch->total_steps-laststridesize) PetscFunctionReturn(0); 671 } 672 if (!stack->solution_only && localstepnum == 0) PetscFunctionReturn(0); /* skip last point in each stride at recompute stage or last stride */ 673 if (stack->solution_only && localstepnum == tjsch->stride-1) PetscFunctionReturn(0); /* skip last step in each stride at recompute stage or last stride */ 674 675 ierr = ElementCreate(ts,stack,&e);CHKERRQ(ierr); 676 ierr = ElementSet(ts,stack,&e,stepnum,time,X);CHKERRQ(ierr); 677 ierr = StackPush(stack,e);CHKERRQ(ierr); 678 PetscFunctionReturn(0); 679 } 680 681 static PetscErrorCode GetTrajTLNR(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum) 682 { 683 Stack *stack = &tjsch->stack; 684 PetscInt id,localstepnum,laststridesize; 685 StackElement e; 686 PetscErrorCode ierr; 687 688 PetscFunctionBegin; 689 if (stepnum == tjsch->total_steps) { 690 ierr = TurnBackward(ts);CHKERRQ(ierr); 691 PetscFunctionReturn(0); 692 } 693 694 localstepnum = stepnum%tjsch->stride; 695 laststridesize = tjsch->total_steps%tjsch->stride; 696 if (!laststridesize) laststridesize = tjsch->stride; 697 if (stack->solution_only) { 698 /* fill stack with info */ 699 if (localstepnum == 0 && tjsch->total_steps-stepnum >= laststridesize) { 700 id = stepnum/tjsch->stride; 701 if (tjsch->save_stack) { 702 ierr = StackLoadAll(tj,ts,stack,id);CHKERRQ(ierr); 703 tjsch->recompute = PETSC_TRUE; 704 tjsch->skip_trajectory = PETSC_TRUE; 705 ierr = TurnForward(ts);CHKERRQ(ierr); 706 ierr = ReCompute(ts,tjsch,id*tjsch->stride-1,id*tjsch->stride);CHKERRQ(ierr); 707 tjsch->skip_trajectory = PETSC_FALSE; 708 } else { 709 ierr = LoadSingle(tj,ts,stack,id);CHKERRQ(ierr); 710 tjsch->recompute = PETSC_TRUE; 711 ierr = TurnForward(ts);CHKERRQ(ierr); 712 ierr = ReCompute(ts,tjsch,(id-1)*tjsch->stride,id*tjsch->stride);CHKERRQ(ierr); 713 } 714 PetscFunctionReturn(0); 715 } 716 /* restore a checkpoint */ 717 ierr = StackPop(stack,&e);CHKERRQ(ierr); 718 ierr = UpdateTS(ts,stack,e);CHKERRQ(ierr); 719 tjsch->recompute = PETSC_TRUE; 720 tjsch->skip_trajectory = PETSC_TRUE; 721 ierr = TurnForward(ts);CHKERRQ(ierr); 722 ierr = ReCompute(ts,tjsch,e->stepnum,stepnum);CHKERRQ(ierr); 723 tjsch->skip_trajectory = PETSC_FALSE; 724 ierr = ElementDestroy(stack,e);CHKERRQ(ierr); 725 } else { 726 /* fill stack with info */ 727 if (localstepnum == 0 && tjsch->total_steps-stepnum >= laststridesize) { 728 id = stepnum/tjsch->stride; 729 if (tjsch->save_stack) { 730 ierr = StackLoadAll(tj,ts,stack,id);CHKERRQ(ierr); 731 } else { 732 ierr = LoadSingle(tj,ts,stack,id);CHKERRQ(ierr); 733 ierr = ElementCreate(ts,stack,&e);CHKERRQ(ierr); 734 ierr = ElementSet(ts,stack,&e,(id-1)*tjsch->stride+1,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 735 ierr = StackPush(stack,e);CHKERRQ(ierr); 736 tjsch->recompute = PETSC_TRUE; 737 ierr = TurnForward(ts);CHKERRQ(ierr); 738 ierr = ReCompute(ts,tjsch,e->stepnum,id*tjsch->stride);CHKERRQ(ierr); 739 } 740 PetscFunctionReturn(0); 741 } 742 /* restore a checkpoint */ 743 ierr = StackPop(stack,&e);CHKERRQ(ierr); 744 ierr = UpdateTS(ts,stack,e);CHKERRQ(ierr); 745 ierr = ElementDestroy(stack,e);CHKERRQ(ierr); 746 } 747 PetscFunctionReturn(0); 748 } 749 750 #if defined(PETSC_HAVE_REVOLVE) 751 static PetscErrorCode printwhattodo(PetscViewer viewer,PetscInt whattodo,RevolveCTX *rctx,PetscInt shift) 752 { 753 PetscErrorCode ierr; 754 755 PetscFunctionBegin; 756 if (!viewer) PetscFunctionReturn(0); 757 758 switch(whattodo) { 759 case 1: 760 ierr = PetscViewerASCIIPrintf(viewer,"Advance from %D to %D\n",rctx->oldcapo+shift,rctx->capo+shift);CHKERRQ(ierr); 761 break; 762 case 2: 763 ierr = PetscViewerASCIIPrintf(viewer,"Store in checkpoint number %D (located in RAM)\n",rctx->check);CHKERRQ(ierr); 764 break; 765 case 3: 766 ierr = PetscViewerASCIIPrintf(viewer,"First turn: Initialize adjoints and reverse first step\n");CHKERRQ(ierr); 767 break; 768 case 4: 769 ierr = PetscViewerASCIIPrintf(viewer,"Forward and reverse one step\n");CHKERRQ(ierr); 770 break; 771 case 5: 772 ierr = PetscViewerASCIIPrintf(viewer,"Restore in checkpoint number %D (located in RAM)\n",rctx->check);CHKERRQ(ierr); 773 break; 774 case 7: 775 ierr = PetscViewerASCIIPrintf(viewer,"Store in checkpoint number %D (located on disk)\n",rctx->check);CHKERRQ(ierr); 776 break; 777 case 8: 778 ierr = PetscViewerASCIIPrintf(viewer,"Restore in checkpoint number %D (located on disk)\n",rctx->check);CHKERRQ(ierr); 779 break; 780 case -1: 781 ierr = PetscViewerASCIIPrintf(viewer,"Error!");CHKERRQ(ierr); 782 break; 783 } 784 PetscFunctionReturn(0); 785 } 786 787 static PetscErrorCode printwhattodo2(PetscViewer viewer,PetscInt whattodo,RevolveCTX *rctx,PetscInt shift) 788 { 789 PetscErrorCode ierr; 790 791 PetscFunctionBegin; 792 if (!viewer) PetscFunctionReturn(0); 793 794 switch(whattodo) { 795 case 1: 796 ierr = PetscViewerASCIIPrintf(viewer,"[Top Level] Advance from stride %D to stride %D\n",rctx->oldcapo+shift,rctx->capo+shift);CHKERRQ(ierr); 797 break; 798 case 2: 799 ierr = PetscViewerASCIIPrintf(viewer,"[Top Level] Store in checkpoint number %D\n",rctx->check);CHKERRQ(ierr); 800 break; 801 case 3: 802 ierr = PetscViewerASCIIPrintf(viewer,"[Top Level] First turn: Initialize adjoints and reverse first stride\n");CHKERRQ(ierr); 803 break; 804 case 4: 805 ierr = PetscViewerASCIIPrintf(viewer,"[Top Level] Forward and reverse one stride\n");CHKERRQ(ierr); 806 break; 807 case 5: 808 ierr = PetscViewerASCIIPrintf(viewer,"[Top Level] Restore in checkpoint number %D\n",rctx->check);CHKERRQ(ierr); 809 break; 810 case 7: 811 ierr = PetscViewerASCIIPrintf(viewer,"[Top Level] Store in top-level checkpoint number %D\n",rctx->check);CHKERRQ(ierr); 812 break; 813 case 8: 814 ierr = PetscViewerASCIIPrintf(viewer,"[Top Level] Restore in top-level checkpoint number %D\n",rctx->check);CHKERRQ(ierr); 815 break; 816 case -1: 817 ierr = PetscViewerASCIIPrintf(viewer,"[Top Level] Error!");CHKERRQ(ierr); 818 break; 819 } 820 PetscFunctionReturn(0); 821 } 822 823 static PetscErrorCode InitRevolve(PetscInt fine,PetscInt snaps,RevolveCTX *rctx) 824 { 825 PetscFunctionBegin; 826 revolve_reset(); 827 revolve_create_offline(fine,snaps); 828 rctx->snaps_in = snaps; 829 rctx->fine = fine; 830 rctx->check = 0; 831 rctx->capo = 0; 832 rctx->reverseonestep = PETSC_FALSE; 833 /* check stepsleft? */ 834 PetscFunctionReturn(0); 835 } 836 837 static PetscErrorCode FastForwardRevolve(RevolveCTX *rctx) 838 { 839 PetscInt whattodo; 840 841 PetscFunctionBegin; 842 whattodo = 0; 843 while(whattodo!=3) { /* we have to fast forward revolve to the beginning of the backward sweep due to unfriendly revolve interface */ 844 whattodo = revolve_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where); 845 } 846 PetscFunctionReturn(0); 847 } 848 849 static PetscErrorCode ApplyRevolve(PetscViewer viewer,SchedulerType stype,RevolveCTX *rctx,PetscInt total_steps,PetscInt stepnum,PetscInt localstepnum,PetscBool toplevel,PetscInt *store) 850 { 851 PetscErrorCode ierr; 852 PetscInt shift,whattodo; 853 854 PetscFunctionBegin; 855 *store = 0; 856 if (rctx->stepsleft > 0) { /* advance the solution without checkpointing anything as Revolve requires */ 857 rctx->stepsleft--; 858 PetscFunctionReturn(0); 859 } 860 /* let Revolve determine what to do next */ 861 shift = stepnum-localstepnum; 862 rctx->oldcapo = rctx->capo; 863 rctx->capo = localstepnum; 864 865 if (!toplevel) whattodo = revolve_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where); 866 else whattodo = revolve2_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where); 867 if (stype == REVOLVE_ONLINE && whattodo == 8) whattodo = 5; 868 if (stype == REVOLVE_ONLINE && whattodo == 7) whattodo = 2; 869 if (!toplevel) {ierr = printwhattodo(viewer,whattodo,rctx,shift);CHKERRQ(ierr);} 870 else {ierr = printwhattodo2(viewer,whattodo,rctx,shift);CHKERRQ(ierr);} 871 if (whattodo == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in the Revolve library"); 872 if (whattodo == 1) { /* advance some time steps */ 873 if (stype == REVOLVE_ONLINE && rctx->capo >= total_steps-1) { 874 revolve_turn(total_steps,&rctx->capo,&rctx->fine); 875 if (!toplevel) {ierr = printwhattodo(viewer,whattodo,rctx,shift);CHKERRQ(ierr);} 876 else {ierr = printwhattodo2(viewer,whattodo,rctx,shift);CHKERRQ(ierr);} 877 } 878 rctx->stepsleft = rctx->capo-rctx->oldcapo-1; 879 } 880 if (whattodo == 3 || whattodo == 4) { /* ready for a reverse step */ 881 rctx->reverseonestep = PETSC_TRUE; 882 } 883 if (whattodo == 5) { /* restore a checkpoint and ask Revolve what to do next */ 884 rctx->oldcapo = rctx->capo; 885 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*/ 886 else whattodo = revolve2_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where); 887 if (!toplevel) {ierr = printwhattodo(viewer,whattodo,rctx,shift);CHKERRQ(ierr);} 888 else {ierr = printwhattodo2(viewer,whattodo,rctx,shift);CHKERRQ(ierr);} 889 if (whattodo == 3 || whattodo == 4) rctx->reverseonestep = PETSC_TRUE; 890 if (whattodo == 1) rctx->stepsleft = rctx->capo-rctx->oldcapo; 891 } 892 if (whattodo == 7) { /* save the checkpoint to disk */ 893 *store = 2; 894 rctx->oldcapo = rctx->capo; 895 whattodo = revolve_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where); /* must return 1 */ 896 ierr = printwhattodo(viewer,whattodo,rctx,shift);CHKERRQ(ierr); 897 rctx->stepsleft = rctx->capo-rctx->oldcapo-1; 898 } 899 if (whattodo == 2) { /* store a checkpoint to RAM and ask Revolve how many time steps to advance next */ 900 *store = 1; 901 rctx->oldcapo = rctx->capo; 902 if (!toplevel) whattodo = revolve_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where); /* must return 1 */ 903 else whattodo = revolve2_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where); 904 if (!toplevel) {ierr = printwhattodo(viewer,whattodo,rctx,shift);CHKERRQ(ierr);} 905 else {ierr = printwhattodo2(viewer,whattodo,rctx,shift);CHKERRQ(ierr);} 906 if (stype == REVOLVE_ONLINE && rctx->capo >= total_steps-1) { 907 revolve_turn(total_steps,&rctx->capo,&rctx->fine); 908 ierr = printwhattodo(viewer,whattodo,rctx,shift);CHKERRQ(ierr); 909 } 910 rctx->stepsleft = rctx->capo-rctx->oldcapo-1; 911 } 912 PetscFunctionReturn(0); 913 } 914 915 static PetscErrorCode SetTrajROF(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X) 916 { 917 Stack *stack = &tjsch->stack; 918 PetscInt store; 919 StackElement e; 920 PetscErrorCode ierr; 921 922 PetscFunctionBegin; 923 if (!stack->solution_only && stepnum == 0) PetscFunctionReturn(0); 924 if (stack->solution_only && stepnum == tjsch->total_steps) PetscFunctionReturn(0); 925 ierr = ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,tjsch->total_steps,stepnum,stepnum,PETSC_FALSE,&store);CHKERRQ(ierr); 926 if (store == 1) { 927 if (stepnum < stack->top) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_MEMC,"Illegal modification of a non-top stack element"); 928 ierr = ElementCreate(ts,stack,&e);CHKERRQ(ierr); 929 ierr = ElementSet(ts,stack,&e,stepnum,time,X);CHKERRQ(ierr); 930 ierr = StackPush(stack,e);CHKERRQ(ierr); 931 } 932 PetscFunctionReturn(0); 933 } 934 935 static PetscErrorCode GetTrajROF(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum) 936 { 937 Stack *stack = &tjsch->stack; 938 PetscInt whattodo,shift,store; 939 StackElement e; 940 PetscErrorCode ierr; 941 942 PetscFunctionBegin; 943 if (stepnum == 0 || stepnum == tjsch->total_steps) { 944 ierr = TurnBackward(ts);CHKERRQ(ierr); 945 tjsch->rctx->reverseonestep = PETSC_FALSE; 946 PetscFunctionReturn(0); 947 } 948 /* restore a checkpoint */ 949 ierr = StackTop(stack,&e);CHKERRQ(ierr); 950 ierr = UpdateTS(ts,stack,e);CHKERRQ(ierr); 951 if (stack->solution_only) { /* start with restoring a checkpoint */ 952 tjsch->rctx->capo = stepnum; 953 tjsch->rctx->oldcapo = tjsch->rctx->capo; 954 shift = 0; 955 whattodo = revolve_action(&tjsch->rctx->check,&tjsch->rctx->capo,&tjsch->rctx->fine,tjsch->rctx->snaps_in,&tjsch->rctx->info,&tjsch->rctx->where); 956 ierr = printwhattodo(tj->monitor,whattodo,tjsch->rctx,shift);CHKERRQ(ierr); 957 } else { /* 2 revolve actions: restore a checkpoint and then advance */ 958 ierr = ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,tjsch->total_steps,stepnum,stepnum,PETSC_FALSE,&store);CHKERRQ(ierr); 959 if (tj->monitor) { 960 ierr = PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 961 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); 962 ierr = PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 963 } 964 if (!tjsch->rctx->reverseonestep && tjsch->rctx->stepsleft > 0) tjsch->rctx->stepsleft--; 965 } 966 if (stack->solution_only || (!stack->solution_only && e->stepnum < stepnum)) { 967 tjsch->recompute = PETSC_TRUE; 968 ierr = TurnForward(ts);CHKERRQ(ierr); 969 ierr = ReCompute(ts,tjsch,e->stepnum,stepnum);CHKERRQ(ierr); 970 } 971 if ((stack->solution_only && e->stepnum+1 == stepnum) || (!stack->solution_only && e->stepnum == stepnum)) { 972 ierr = StackPop(stack,&e);CHKERRQ(ierr); 973 ierr = ElementDestroy(stack,e);CHKERRQ(ierr); 974 } 975 tjsch->rctx->reverseonestep = PETSC_FALSE; 976 PetscFunctionReturn(0); 977 } 978 979 static PetscErrorCode SetTrajRON(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X) 980 { 981 Stack *stack = &tjsch->stack; 982 Vec *Y; 983 PetscInt i,store; 984 PetscReal timeprev; 985 StackElement e; 986 RevolveCTX *rctx = tjsch->rctx; 987 PetscErrorCode ierr; 988 989 PetscFunctionBegin; 990 if (!stack->solution_only && stepnum == 0) PetscFunctionReturn(0); 991 if (stack->solution_only && stepnum == tjsch->total_steps) PetscFunctionReturn(0); 992 ierr = ApplyRevolve(tj->monitor,tjsch->stype,rctx,tjsch->total_steps,stepnum,stepnum,PETSC_FALSE,&store);CHKERRQ(ierr); 993 if (store == 1) { 994 if (rctx->check != stack->top+1) { /* overwrite some non-top checkpoint in the stack */ 995 ierr = StackFind(stack,&e,rctx->check);CHKERRQ(ierr); 996 ierr = VecCopy(X,e->X);CHKERRQ(ierr); 997 if (stack->numY > 0 && !stack->solution_only) { 998 ierr = TSGetStages(ts,&stack->numY,&Y);CHKERRQ(ierr); 999 for (i=0;i<stack->numY;i++) { 1000 ierr = VecCopy(Y[i],e->Y[i]);CHKERRQ(ierr); 1001 } 1002 } 1003 e->stepnum = stepnum; 1004 e->time = time; 1005 ierr = TSGetPrevTime(ts,&timeprev);CHKERRQ(ierr); 1006 e->timeprev = timeprev; 1007 } else { 1008 if (stepnum < stack->top) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_MEMC,"Illegal modification of a non-top stack element"); 1009 ierr = ElementCreate(ts,stack,&e);CHKERRQ(ierr); 1010 ierr = ElementSet(ts,stack,&e,stepnum,time,X);CHKERRQ(ierr); 1011 ierr = StackPush(stack,e);CHKERRQ(ierr); 1012 } 1013 } 1014 PetscFunctionReturn(0); 1015 } 1016 1017 static PetscErrorCode GetTrajRON(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum) 1018 { 1019 Stack *stack = &tjsch->stack; 1020 PetscInt whattodo,shift; 1021 StackElement e; 1022 PetscErrorCode ierr; 1023 1024 PetscFunctionBegin; 1025 if (stepnum == 0 || stepnum == tjsch->total_steps) { 1026 ierr = TurnBackward(ts);CHKERRQ(ierr); 1027 tjsch->rctx->reverseonestep = PETSC_FALSE; 1028 PetscFunctionReturn(0); 1029 } 1030 tjsch->rctx->capo = stepnum; 1031 tjsch->rctx->oldcapo = tjsch->rctx->capo; 1032 shift = 0; 1033 whattodo = revolve_action(&tjsch->rctx->check,&tjsch->rctx->capo,&tjsch->rctx->fine,tjsch->rctx->snaps_in,&tjsch->rctx->info,&tjsch->rctx->where); /* whattodo=restore */ 1034 if (whattodo == 8) whattodo = 5; 1035 ierr = printwhattodo(tj->monitor,whattodo,tjsch->rctx,shift);CHKERRQ(ierr); 1036 /* restore a checkpoint */ 1037 ierr = StackFind(stack,&e,tjsch->rctx->check);CHKERRQ(ierr); 1038 ierr = UpdateTS(ts,stack,e);CHKERRQ(ierr); 1039 if (!stack->solution_only) { /* whattodo must be 5 */ 1040 /* ask Revolve what to do next */ 1041 tjsch->rctx->oldcapo = tjsch->rctx->capo; 1042 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*/ 1043 ierr = printwhattodo(tj->monitor,whattodo,tjsch->rctx,shift);CHKERRQ(ierr); 1044 if (whattodo == 3 || whattodo == 4) tjsch->rctx->reverseonestep = PETSC_TRUE; 1045 if (whattodo == 1) tjsch->rctx->stepsleft = tjsch->rctx->capo-tjsch->rctx->oldcapo; 1046 if (tj->monitor) { 1047 ierr = PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 1048 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); 1049 ierr = PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 1050 } 1051 if (!tjsch->rctx->reverseonestep && tjsch->rctx->stepsleft > 0) tjsch->rctx->stepsleft--; 1052 } 1053 if (stack->solution_only || (!stack->solution_only && e->stepnum < stepnum)) { 1054 tjsch->recompute = PETSC_TRUE; 1055 ierr = TurnForward(ts);CHKERRQ(ierr); 1056 ierr = ReCompute(ts,tjsch,e->stepnum,stepnum);CHKERRQ(ierr); 1057 } 1058 tjsch->rctx->reverseonestep = PETSC_FALSE; 1059 PetscFunctionReturn(0); 1060 } 1061 1062 static PetscErrorCode SetTrajTLR(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X) 1063 { 1064 Stack *stack = &tjsch->stack; 1065 PetscInt store,localstepnum,laststridesize; 1066 StackElement e; 1067 PetscBool done = PETSC_FALSE; 1068 PetscErrorCode ierr; 1069 1070 PetscFunctionBegin; 1071 if (!stack->solution_only && stepnum == 0) PetscFunctionReturn(0); 1072 if (stack->solution_only && stepnum == tjsch->total_steps) PetscFunctionReturn(0); 1073 1074 localstepnum = stepnum%tjsch->stride; 1075 laststridesize = tjsch->total_steps%tjsch->stride; 1076 if (!laststridesize) laststridesize = tjsch->stride; 1077 1078 if (!tjsch->recompute) { 1079 ierr = TopLevelStore(tj,ts,tjsch,stepnum,localstepnum,laststridesize,&done);CHKERRQ(ierr); 1080 /* revolve is needed for the last stride; different starting points for last stride between solutin_only and !solutin_only */ 1081 if (!stack->solution_only && !tjsch->save_stack && stepnum <= tjsch->total_steps-laststridesize) PetscFunctionReturn(0); 1082 if (stack->solution_only && !tjsch->save_stack && stepnum < tjsch->total_steps-laststridesize) PetscFunctionReturn(0); 1083 } 1084 if (tjsch->save_stack && done) { 1085 ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr); 1086 PetscFunctionReturn(0); 1087 } 1088 if (laststridesize < tjsch->stride) { 1089 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 */ 1090 ierr = InitRevolve(laststridesize,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr); 1091 } 1092 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 */ 1093 ierr = InitRevolve(laststridesize,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr); 1094 } 1095 } 1096 ierr = ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,tjsch->total_steps,stepnum,localstepnum,PETSC_FALSE,&store);CHKERRQ(ierr); 1097 if (store == 1) { 1098 if (localstepnum < stack->top) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_MEMC,"Illegal modification of a non-top stack element"); 1099 ierr = ElementCreate(ts,stack,&e);CHKERRQ(ierr); 1100 ierr = ElementSet(ts,stack,&e,stepnum,time,X);CHKERRQ(ierr); 1101 ierr = StackPush(stack,e);CHKERRQ(ierr); 1102 } 1103 PetscFunctionReturn(0); 1104 } 1105 1106 static PetscErrorCode GetTrajTLR(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum) 1107 { 1108 Stack *stack = &tjsch->stack; 1109 PetscInt whattodo,shift; 1110 PetscInt localstepnum,stridenum,laststridesize,store; 1111 StackElement e; 1112 PetscErrorCode ierr; 1113 1114 PetscFunctionBegin; 1115 localstepnum = stepnum%tjsch->stride; 1116 stridenum = stepnum/tjsch->stride; 1117 if (stepnum == tjsch->total_steps) { 1118 ierr = TurnBackward(ts);CHKERRQ(ierr); 1119 tjsch->rctx->reverseonestep = PETSC_FALSE; 1120 PetscFunctionReturn(0); 1121 } 1122 laststridesize = tjsch->total_steps%tjsch->stride; 1123 if (!laststridesize) laststridesize = tjsch->stride; 1124 if (stack->solution_only) { 1125 /* fill stack */ 1126 if (localstepnum == 0 && stepnum <= tjsch->total_steps-laststridesize) { 1127 if (tjsch->save_stack) { 1128 ierr = StackLoadAll(tj,ts,stack,stridenum);CHKERRQ(ierr); 1129 ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr); 1130 ierr = FastForwardRevolve(tjsch->rctx);CHKERRQ(ierr); 1131 tjsch->recompute = PETSC_TRUE; 1132 tjsch->skip_trajectory = PETSC_TRUE; 1133 ierr = TurnForward(ts);CHKERRQ(ierr); 1134 ierr = ReCompute(ts,tjsch,stridenum*tjsch->stride-1,stridenum*tjsch->stride);CHKERRQ(ierr); 1135 tjsch->skip_trajectory = PETSC_FALSE; 1136 } else { 1137 ierr = LoadSingle(tj,ts,stack,stridenum);CHKERRQ(ierr); 1138 ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr); 1139 tjsch->recompute = PETSC_TRUE; 1140 ierr = TurnForward(ts);CHKERRQ(ierr); 1141 ierr = ReCompute(ts,tjsch,(stridenum-1)*tjsch->stride,stridenum*tjsch->stride);CHKERRQ(ierr); 1142 } 1143 PetscFunctionReturn(0); 1144 } 1145 /* restore a checkpoint */ 1146 ierr = StackTop(stack,&e);CHKERRQ(ierr); 1147 ierr = UpdateTS(ts,stack,e);CHKERRQ(ierr); 1148 /* start with restoring a checkpoint */ 1149 tjsch->rctx->capo = stepnum; 1150 tjsch->rctx->oldcapo = tjsch->rctx->capo; 1151 shift = stepnum-localstepnum; 1152 whattodo = revolve_action(&tjsch->rctx->check,&tjsch->rctx->capo,&tjsch->rctx->fine,tjsch->rctx->snaps_in,&tjsch->rctx->info,&tjsch->rctx->where); 1153 ierr = printwhattodo(tj->monitor,whattodo,tjsch->rctx,shift);CHKERRQ(ierr); 1154 tjsch->recompute = PETSC_TRUE; 1155 ierr = TurnForward(ts);CHKERRQ(ierr); 1156 ierr = ReCompute(ts,tjsch,e->stepnum,stepnum);CHKERRQ(ierr); 1157 if (e->stepnum+1 == stepnum) { 1158 ierr = StackPop(stack,&e);CHKERRQ(ierr); 1159 ierr = ElementDestroy(stack,e);CHKERRQ(ierr); 1160 } 1161 } else { 1162 /* fill stack with info */ 1163 if (localstepnum == 0 && tjsch->total_steps-stepnum >= laststridesize) { 1164 if (tjsch->save_stack) { 1165 ierr = StackLoadAll(tj,ts,stack,stridenum);CHKERRQ(ierr); 1166 ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr); 1167 ierr = FastForwardRevolve(tjsch->rctx);CHKERRQ(ierr); 1168 } else { 1169 ierr = LoadSingle(tj,ts,stack,stridenum);CHKERRQ(ierr); 1170 ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr); 1171 ierr = ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,tjsch->total_steps,(stridenum-1)*tjsch->stride+1,1,PETSC_FALSE,&store);CHKERRQ(ierr); 1172 if (tj->monitor) { 1173 ierr = PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 1174 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); 1175 ierr = PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 1176 } 1177 ierr = ElementCreate(ts,stack,&e);CHKERRQ(ierr); 1178 ierr = ElementSet(ts,stack,&e,(stridenum-1)*tjsch->stride+1,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 1179 ierr = StackPush(stack,e);CHKERRQ(ierr); 1180 tjsch->recompute = PETSC_TRUE; 1181 ierr = TurnForward(ts);CHKERRQ(ierr); 1182 ierr = ReCompute(ts,tjsch,e->stepnum,stridenum*tjsch->stride);CHKERRQ(ierr); 1183 } 1184 PetscFunctionReturn(0); 1185 } 1186 /* restore a checkpoint */ 1187 ierr = StackTop(stack,&e);CHKERRQ(ierr); 1188 ierr = UpdateTS(ts,stack,e);CHKERRQ(ierr); 1189 /* 2 revolve actions: restore a checkpoint and then advance */ 1190 ierr = ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,tjsch->total_steps,stepnum,localstepnum,PETSC_FALSE,&store);CHKERRQ(ierr); 1191 if (tj->monitor) { 1192 ierr = PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 1193 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); 1194 ierr = PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 1195 } 1196 if (!tjsch->rctx->reverseonestep && tjsch->rctx->stepsleft > 0) tjsch->rctx->stepsleft--; 1197 if (e->stepnum < stepnum) { 1198 tjsch->recompute = PETSC_TRUE; 1199 ierr = TurnForward(ts);CHKERRQ(ierr); 1200 ierr = ReCompute(ts,tjsch,e->stepnum,stepnum);CHKERRQ(ierr); 1201 } 1202 if (e->stepnum == stepnum) { 1203 ierr = StackPop(stack,&e);CHKERRQ(ierr); 1204 ierr = ElementDestroy(stack,e);CHKERRQ(ierr); 1205 } 1206 } 1207 tjsch->rctx->reverseonestep = PETSC_FALSE; 1208 PetscFunctionReturn(0); 1209 } 1210 1211 static PetscErrorCode SetTrajTLTR(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X) 1212 { 1213 Stack *stack = &tjsch->stack; 1214 PetscInt store,localstepnum,stridenum,laststridesize; 1215 StackElement e; 1216 PetscBool done = PETSC_FALSE; 1217 PetscErrorCode ierr; 1218 1219 PetscFunctionBegin; 1220 if (!stack->solution_only && stepnum == 0) PetscFunctionReturn(0); 1221 if (stack->solution_only && stepnum == tjsch->total_steps) PetscFunctionReturn(0); 1222 1223 localstepnum = stepnum%tjsch->stride; /* index at the bottom level (inside a stride) */ 1224 stridenum = stepnum/tjsch->stride; /* index at the top level */ 1225 laststridesize = tjsch->total_steps%tjsch->stride; 1226 if (!laststridesize) laststridesize = tjsch->stride; 1227 if (stack->solution_only && localstepnum == 0 && !tjsch->rctx2->reverseonestep) { 1228 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); 1229 if (laststridesize < tjsch->stride && stepnum == tjsch->total_steps-laststridesize) { 1230 ierr = InitRevolve(laststridesize,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr); 1231 } 1232 } 1233 if (!stack->solution_only && localstepnum == 1 && !tjsch->rctx2->reverseonestep) { 1234 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); 1235 if (laststridesize < tjsch->stride && stepnum == tjsch->total_steps-laststridesize+1) { 1236 ierr = InitRevolve(laststridesize,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr); 1237 } 1238 } 1239 if (tjsch->store_stride) { 1240 ierr = TopLevelStore(tj,ts,tjsch,stepnum,localstepnum,laststridesize,&done);CHKERRQ(ierr); 1241 if (done) { 1242 ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr); 1243 PetscFunctionReturn(0); 1244 } 1245 } 1246 if (stepnum < tjsch->total_steps-laststridesize) { 1247 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 */ 1248 if (!tjsch->save_stack && !tjsch->rctx2->reverseonestep) PetscFunctionReturn(0); /* store operation does not require revolve be called at bottom level */ 1249 } 1250 /* 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 */ 1251 if (!stack->solution_only && localstepnum == 0 && stepnum != tjsch->total_steps && !tjsch->recompute) PetscFunctionReturn(0); 1252 ierr = ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,tjsch->total_steps,stepnum,localstepnum,PETSC_FALSE,&store);CHKERRQ(ierr); 1253 if (store == 1) { 1254 if (localstepnum < stack->top) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_MEMC,"Illegal modification of a non-top stack element"); 1255 ierr = ElementCreate(ts,stack,&e);CHKERRQ(ierr); 1256 ierr = ElementSet(ts,stack,&e,stepnum,time,X);CHKERRQ(ierr); 1257 ierr = StackPush(stack,e);CHKERRQ(ierr); 1258 } 1259 PetscFunctionReturn(0); 1260 } 1261 1262 static PetscErrorCode GetTrajTLTR(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum) 1263 { 1264 Stack *stack = &tjsch->stack; 1265 DiskStack *diskstack = &tjsch->diskstack; 1266 PetscInt whattodo,shift; 1267 PetscInt localstepnum,stridenum,restoredstridenum,laststridesize,store; 1268 StackElement e; 1269 PetscErrorCode ierr; 1270 1271 PetscFunctionBegin; 1272 localstepnum = stepnum%tjsch->stride; 1273 stridenum = stepnum/tjsch->stride; 1274 if (stepnum == tjsch->total_steps) { 1275 ierr = TurnBackward(ts);CHKERRQ(ierr); 1276 tjsch->rctx->reverseonestep = PETSC_FALSE; 1277 PetscFunctionReturn(0); 1278 } 1279 laststridesize = tjsch->total_steps%tjsch->stride; 1280 if (!laststridesize) laststridesize = tjsch->stride; 1281 /* 1282 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: 1283 Case 1 (save_stack) 1284 Restore a disk checkpoint; update TS with the last element in the restored data; recompute to the current point. 1285 Case 2 (!save_stack) 1286 Restore a disk checkpoint; update TS with the restored point; recompute to the current point. 1287 */ 1288 if (localstepnum == 0 && stepnum <= tjsch->total_steps-laststridesize) { 1289 /* restore the top element in the stack for disk checkpoints */ 1290 restoredstridenum = diskstack->container[diskstack->top]; 1291 tjsch->rctx2->reverseonestep = PETSC_FALSE; 1292 /* top-level revolve must be applied before current step, just like the solution_only mode for single-level revolve */ 1293 if (!tjsch->save_stack && stack->solution_only) { /* start with restoring a checkpoint */ 1294 tjsch->rctx2->capo = stridenum; 1295 tjsch->rctx2->oldcapo = tjsch->rctx2->capo; 1296 shift = 0; 1297 whattodo = revolve2_action(&tjsch->rctx2->check,&tjsch->rctx2->capo,&tjsch->rctx2->fine,tjsch->rctx2->snaps_in,&tjsch->rctx2->info,&tjsch->rctx2->where); 1298 ierr = printwhattodo2(tj->monitor,whattodo,tjsch->rctx2,shift);CHKERRQ(ierr); 1299 } else { /* 2 revolve actions: restore a checkpoint and then advance */ 1300 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); 1301 if (tj->monitor) { 1302 ierr = PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 1303 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); 1304 ierr = PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 1305 } 1306 if (!tjsch->rctx2->reverseonestep && tjsch->rctx2->stepsleft > 0) tjsch->rctx2->stepsleft--; 1307 } 1308 /* fill stack */ 1309 if (stack->solution_only) { 1310 if (tjsch->save_stack) { 1311 if (restoredstridenum < stridenum) { 1312 ierr = StackLoadLast(tj,ts,stack,restoredstridenum);CHKERRQ(ierr); 1313 } else { 1314 ierr = StackLoadAll(tj,ts,stack,restoredstridenum);CHKERRQ(ierr); 1315 } 1316 /* recompute one step ahead */ 1317 tjsch->recompute = PETSC_TRUE; 1318 tjsch->skip_trajectory = PETSC_TRUE; 1319 ierr = TurnForward(ts);CHKERRQ(ierr); 1320 ierr = ReCompute(ts,tjsch,stridenum*tjsch->stride-1,stridenum*tjsch->stride);CHKERRQ(ierr); 1321 tjsch->skip_trajectory = PETSC_FALSE; 1322 if (restoredstridenum < stridenum) { 1323 ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr); 1324 tjsch->recompute = PETSC_TRUE; 1325 ierr = TurnForward(ts);CHKERRQ(ierr); 1326 ierr = ReCompute(ts,tjsch,restoredstridenum*tjsch->stride,stepnum);CHKERRQ(ierr); 1327 } else { /* stack ready, fast forward revolve status */ 1328 ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr); 1329 ierr = FastForwardRevolve(tjsch->rctx);CHKERRQ(ierr); 1330 } 1331 } else { 1332 ierr = LoadSingle(tj,ts,stack,restoredstridenum);CHKERRQ(ierr); 1333 ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr); 1334 tjsch->recompute = PETSC_TRUE; 1335 ierr = TurnForward(ts);CHKERRQ(ierr); 1336 ierr = ReCompute(ts,tjsch,(restoredstridenum-1)*tjsch->stride,stepnum);CHKERRQ(ierr); 1337 } 1338 } else { 1339 if (tjsch->save_stack) { 1340 if (restoredstridenum < stridenum) { 1341 ierr = StackLoadLast(tj,ts,stack,restoredstridenum);CHKERRQ(ierr); 1342 /* reset revolve */ 1343 ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr); 1344 tjsch->recompute = PETSC_TRUE; 1345 ierr = TurnForward(ts);CHKERRQ(ierr); 1346 ierr = ReCompute(ts,tjsch,restoredstridenum*tjsch->stride,stepnum);CHKERRQ(ierr); 1347 } else { /* stack ready, fast forward revolve status */ 1348 ierr = StackLoadAll(tj,ts,stack,restoredstridenum);CHKERRQ(ierr); 1349 ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr); 1350 ierr = FastForwardRevolve(tjsch->rctx);CHKERRQ(ierr); 1351 } 1352 } else { 1353 ierr = LoadSingle(tj,ts,stack,restoredstridenum);CHKERRQ(ierr); 1354 ierr = InitRevolve(tjsch->stride,tjsch->max_cps_ram,tjsch->rctx);CHKERRQ(ierr); 1355 /* push first element to stack */ 1356 if (tjsch->store_stride || tjsch->rctx2->reverseonestep) { 1357 shift = (restoredstridenum-1)*tjsch->stride-localstepnum; 1358 ierr = ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,tjsch->total_steps,(restoredstridenum-1)*tjsch->stride+1,1,PETSC_FALSE,&store);CHKERRQ(ierr); 1359 if (tj->monitor) { 1360 ierr = PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 1361 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); 1362 ierr = PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 1363 } 1364 ierr = ElementCreate(ts,stack,&e);CHKERRQ(ierr); 1365 ierr = ElementSet(ts,stack,&e,(restoredstridenum-1)*tjsch->stride+1,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 1366 ierr = StackPush(stack,e);CHKERRQ(ierr); 1367 } 1368 tjsch->recompute = PETSC_TRUE; 1369 ierr = TurnForward(ts);CHKERRQ(ierr); 1370 ierr = ReCompute(ts,tjsch,(restoredstridenum-1)*tjsch->stride+1,stepnum);CHKERRQ(ierr); 1371 } 1372 } 1373 if (restoredstridenum == stridenum) diskstack->top--; 1374 tjsch->rctx->reverseonestep = PETSC_FALSE; 1375 PetscFunctionReturn(0); 1376 } 1377 1378 if (stack->solution_only) { 1379 /* restore a checkpoint */ 1380 ierr = StackTop(stack,&e);CHKERRQ(ierr); 1381 ierr = UpdateTS(ts,stack,e);CHKERRQ(ierr); 1382 /* start with restoring a checkpoint */ 1383 tjsch->rctx->capo = stepnum; 1384 tjsch->rctx->oldcapo = tjsch->rctx->capo; 1385 shift = stepnum-localstepnum; 1386 whattodo = revolve_action(&tjsch->rctx->check,&tjsch->rctx->capo,&tjsch->rctx->fine,tjsch->rctx->snaps_in,&tjsch->rctx->info,&tjsch->rctx->where); 1387 ierr = printwhattodo(tj->monitor,whattodo,tjsch->rctx,shift);CHKERRQ(ierr); 1388 tjsch->recompute = PETSC_TRUE; 1389 ierr = TurnForward(ts);CHKERRQ(ierr); 1390 ierr = ReCompute(ts,tjsch,e->stepnum,stepnum);CHKERRQ(ierr); 1391 if (e->stepnum+1 == stepnum) { 1392 ierr = StackPop(stack,&e);CHKERRQ(ierr); 1393 ierr = ElementDestroy(stack,e);CHKERRQ(ierr); 1394 } 1395 } else { 1396 /* restore a checkpoint */ 1397 ierr = StackTop(stack,&e);CHKERRQ(ierr); 1398 ierr = UpdateTS(ts,stack,e);CHKERRQ(ierr); 1399 /* 2 revolve actions: restore a checkpoint and then advance */ 1400 ierr = ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,tjsch->total_steps,stepnum,localstepnum,PETSC_FALSE,&store);CHKERRQ(ierr); 1401 if (tj->monitor) { 1402 ierr = PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 1403 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); 1404 ierr = PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 1405 } 1406 if (!tjsch->rctx->reverseonestep && tjsch->rctx->stepsleft > 0) tjsch->rctx->stepsleft--; 1407 if (e->stepnum < stepnum) { 1408 tjsch->recompute = PETSC_TRUE; 1409 ierr = TurnForward(ts);CHKERRQ(ierr); 1410 ierr = ReCompute(ts,tjsch,e->stepnum,stepnum);CHKERRQ(ierr); 1411 } 1412 if (e->stepnum == stepnum) { 1413 ierr = StackPop(stack,&e);CHKERRQ(ierr); 1414 ierr = ElementDestroy(stack,e);CHKERRQ(ierr); 1415 } 1416 } 1417 tjsch->rctx->reverseonestep = PETSC_FALSE; 1418 PetscFunctionReturn(0); 1419 } 1420 1421 static PetscErrorCode SetTrajRMS(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X) 1422 { 1423 Stack *stack = &tjsch->stack; 1424 PetscInt store; 1425 StackElement e; 1426 PetscErrorCode ierr; 1427 1428 PetscFunctionBegin; 1429 if (!stack->solution_only && stepnum == 0) PetscFunctionReturn(0); 1430 if (stack->solution_only && stepnum == tjsch->total_steps) PetscFunctionReturn(0); 1431 ierr = ApplyRevolve(tj->monitor,tjsch->stype,tjsch->rctx,tjsch->total_steps,stepnum,stepnum,PETSC_FALSE,&store);CHKERRQ(ierr); 1432 if (store == 1){ 1433 if (stepnum < stack->top) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_MEMC,"Illegal modification of a non-top stack element"); 1434 ierr = ElementCreate(ts,stack,&e);CHKERRQ(ierr); 1435 ierr = ElementSet(ts,stack,&e,stepnum,time,X);CHKERRQ(ierr); 1436 ierr = StackPush(stack,e);CHKERRQ(ierr); 1437 } else if (store == 2) { 1438 ierr = DumpSingle(tj,ts,stack,tjsch->rctx->check+1);CHKERRQ(ierr); 1439 } 1440 PetscFunctionReturn(0); 1441 } 1442 1443 static PetscErrorCode GetTrajRMS(TSTrajectory tj,TS ts,TJScheduler *tjsch,PetscInt stepnum) 1444 { 1445 Stack *stack = &tjsch->stack; 1446 PetscInt whattodo,shift; 1447 PetscInt restart; 1448 PetscBool ondisk; 1449 StackElement e; 1450 PetscErrorCode ierr; 1451 1452 PetscFunctionBegin; 1453 if (stepnum == 0 || stepnum == tjsch->total_steps) { 1454 ierr = TurnBackward(ts);CHKERRQ(ierr); 1455 tjsch->rctx->reverseonestep = PETSC_FALSE; 1456 PetscFunctionReturn(0); 1457 } 1458 tjsch->rctx->capo = stepnum; 1459 tjsch->rctx->oldcapo = tjsch->rctx->capo; 1460 shift = 0; 1461 whattodo = revolve_action(&tjsch->rctx->check,&tjsch->rctx->capo,&tjsch->rctx->fine,tjsch->rctx->snaps_in,&tjsch->rctx->info,&tjsch->rctx->where); /* whattodo=restore */ 1462 ierr = printwhattodo(tj->monitor,whattodo,tjsch->rctx,shift);CHKERRQ(ierr); 1463 /* restore a checkpoint */ 1464 restart = tjsch->rctx->capo; 1465 if (!tjsch->rctx->where) { 1466 ondisk = PETSC_TRUE; 1467 ierr = LoadSingle(tj,ts,stack,tjsch->rctx->check+1);CHKERRQ(ierr); 1468 ierr = TurnBackward(ts);CHKERRQ(ierr); 1469 } else { 1470 ondisk = PETSC_FALSE; 1471 ierr = StackTop(stack,&e);CHKERRQ(ierr); 1472 ierr = UpdateTS(ts,stack,e);CHKERRQ(ierr); 1473 } 1474 if (!stack->solution_only) { /* whattodo must be 5 or 8 */ 1475 /* ask Revolve what to do next */ 1476 tjsch->rctx->oldcapo = tjsch->rctx->capo; 1477 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*/ 1478 ierr = printwhattodo(tj->monitor,whattodo,tjsch->rctx,shift);CHKERRQ(ierr); 1479 if (whattodo == 3 || whattodo == 4) tjsch->rctx->reverseonestep = PETSC_TRUE; 1480 if (whattodo == 1) tjsch->rctx->stepsleft = tjsch->rctx->capo-tjsch->rctx->oldcapo; 1481 if (tj->monitor) { 1482 ierr = PetscViewerASCIIAddTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 1483 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); 1484 ierr = PetscViewerASCIISubtractTab(tj->monitor,((PetscObject)tj)->tablevel);CHKERRQ(ierr); 1485 } 1486 if (!tjsch->rctx->reverseonestep && tjsch->rctx->stepsleft > 0) tjsch->rctx->stepsleft--; 1487 restart++; /* skip one step */ 1488 } 1489 if (stack->solution_only || (!stack->solution_only && restart < stepnum)) { 1490 tjsch->recompute = PETSC_TRUE; 1491 ierr = TurnForward(ts);CHKERRQ(ierr); 1492 ierr = ReCompute(ts,tjsch,restart,stepnum);CHKERRQ(ierr); 1493 } 1494 if (!ondisk && ( (stack->solution_only && e->stepnum+1 == stepnum) || (!stack->solution_only && e->stepnum == stepnum) )) { 1495 ierr = StackPop(stack,&e);CHKERRQ(ierr); 1496 ierr = ElementDestroy(stack,e);CHKERRQ(ierr); 1497 } 1498 tjsch->rctx->reverseonestep = PETSC_FALSE; 1499 PetscFunctionReturn(0); 1500 } 1501 #endif 1502 1503 static PetscErrorCode TSTrajectorySet_Memory(TSTrajectory tj,TS ts,PetscInt stepnum,PetscReal time,Vec X) 1504 { 1505 TJScheduler *tjsch = (TJScheduler*)tj->data; 1506 PetscErrorCode ierr; 1507 1508 PetscFunctionBegin; 1509 if (!tjsch->recompute) { /* use global stepnum in the forward sweep */ 1510 ierr = TSGetStepNumber(ts,&stepnum);CHKERRQ(ierr); 1511 } 1512 /* for consistency */ 1513 if (!tjsch->recompute && stepnum == 0) ts->ptime_prev = ts->ptime-ts->time_step; 1514 switch (tjsch->stype) { 1515 case NONE: 1516 ierr = SetTrajN(ts,tjsch,stepnum,time,X);CHKERRQ(ierr); 1517 break; 1518 case TWO_LEVEL_NOREVOLVE: 1519 ierr = SetTrajTLNR(tj,ts,tjsch,stepnum,time,X);CHKERRQ(ierr); 1520 break; 1521 #if defined(PETSC_HAVE_REVOLVE) 1522 case TWO_LEVEL_REVOLVE: 1523 ierr = SetTrajTLR(tj,ts,tjsch,stepnum,time,X);CHKERRQ(ierr); 1524 break; 1525 case TWO_LEVEL_TWO_REVOLVE: 1526 ierr = SetTrajTLTR(tj,ts,tjsch,stepnum,time,X);CHKERRQ(ierr); 1527 break; 1528 case REVOLVE_OFFLINE: 1529 ierr = SetTrajROF(tj,ts,tjsch,stepnum,time,X);CHKERRQ(ierr); 1530 break; 1531 case REVOLVE_ONLINE: 1532 ierr = SetTrajRON(tj,ts,tjsch,stepnum,time,X);CHKERRQ(ierr); 1533 break; 1534 case REVOLVE_MULTISTAGE: 1535 ierr = SetTrajRMS(tj,ts,tjsch,stepnum,time,X);CHKERRQ(ierr); 1536 break; 1537 #endif 1538 default: 1539 break; 1540 } 1541 PetscFunctionReturn(0); 1542 } 1543 1544 static PetscErrorCode TSTrajectoryGet_Memory(TSTrajectory tj,TS ts,PetscInt stepnum,PetscReal *t) 1545 { 1546 TJScheduler *tjsch = (TJScheduler*)tj->data; 1547 PetscErrorCode ierr; 1548 1549 PetscFunctionBegin; 1550 if (stepnum == 0) PetscFunctionReturn(0); 1551 switch (tjsch->stype) { 1552 case NONE: 1553 ierr = GetTrajN(ts,tjsch,stepnum);CHKERRQ(ierr); 1554 break; 1555 case TWO_LEVEL_NOREVOLVE: 1556 ierr = GetTrajTLNR(tj,ts,tjsch,stepnum);CHKERRQ(ierr); 1557 break; 1558 #if defined(PETSC_HAVE_REVOLVE) 1559 case TWO_LEVEL_REVOLVE: 1560 ierr = GetTrajTLR(tj,ts,tjsch,stepnum);CHKERRQ(ierr); 1561 break; 1562 case TWO_LEVEL_TWO_REVOLVE: 1563 ierr = GetTrajTLTR(tj,ts,tjsch,stepnum);CHKERRQ(ierr); 1564 break; 1565 case REVOLVE_OFFLINE: 1566 ierr = GetTrajROF(tj,ts,tjsch,stepnum);CHKERRQ(ierr); 1567 break; 1568 case REVOLVE_ONLINE: 1569 ierr = GetTrajRON(tj,ts,tjsch,stepnum);CHKERRQ(ierr); 1570 break; 1571 case REVOLVE_MULTISTAGE: 1572 ierr = GetTrajRMS(tj,ts,tjsch,stepnum);CHKERRQ(ierr); 1573 break; 1574 #endif 1575 default: 1576 break; 1577 } 1578 PetscFunctionReturn(0); 1579 } 1580 1581 PETSC_UNUSED static PetscErrorCode TSTrajectorySetStride_Memory(TSTrajectory tj,PetscInt stride) 1582 { 1583 TJScheduler *tjsch = (TJScheduler*)tj->data; 1584 1585 PetscFunctionBegin; 1586 tjsch->stride = stride; 1587 PetscFunctionReturn(0); 1588 } 1589 1590 PETSC_UNUSED static PetscErrorCode TSTrajectorySetMaxCpsRAM_Memory(TSTrajectory tj,PetscInt max_cps_ram) 1591 { 1592 TJScheduler *tjsch = (TJScheduler*)tj->data; 1593 1594 PetscFunctionBegin; 1595 tjsch->max_cps_ram = max_cps_ram; 1596 PetscFunctionReturn(0); 1597 } 1598 1599 PETSC_UNUSED static PetscErrorCode TSTrajectorySetMaxCpsDisk_Memory(TSTrajectory tj,PetscInt max_cps_disk) 1600 { 1601 TJScheduler *tjsch = (TJScheduler*)tj->data; 1602 1603 PetscFunctionBegin; 1604 tjsch->max_cps_disk = max_cps_disk; 1605 PetscFunctionReturn(0); 1606 } 1607 1608 #if defined(PETSC_HAVE_REVOLVE) 1609 PETSC_UNUSED static PetscErrorCode TSTrajectorySetRevolveOnline(TSTrajectory tj,PetscBool use_online) 1610 { 1611 TJScheduler *tjsch = (TJScheduler*)tj->data; 1612 1613 PetscFunctionBegin; 1614 tjsch->use_online = use_online; 1615 PetscFunctionReturn(0); 1616 } 1617 #endif 1618 1619 PETSC_UNUSED static PetscErrorCode TSTrajectorySetSaveStack(TSTrajectory tj,PetscBool save_stack) 1620 { 1621 TJScheduler *tjsch = (TJScheduler*)tj->data; 1622 1623 PetscFunctionBegin; 1624 tjsch->save_stack = save_stack; 1625 PetscFunctionReturn(0); 1626 } 1627 1628 PETSC_UNUSED static PetscErrorCode TSTrajectorySetSolutionOnly(TSTrajectory tj,PetscBool solution_only) 1629 { 1630 TJScheduler *tjsch = (TJScheduler*)tj->data; 1631 Stack *stack = &tjsch->stack; 1632 1633 PetscFunctionBegin; 1634 stack->solution_only = solution_only; 1635 PetscFunctionReturn(0); 1636 } 1637 1638 PETSC_UNUSED static PetscErrorCode TSTrajectorySetUseDRAM(TSTrajectory tj,PetscBool use_dram) 1639 { 1640 TJScheduler *tjsch = (TJScheduler*)tj->data; 1641 1642 PetscFunctionBegin; 1643 tjsch->stack.use_dram = use_dram; 1644 PetscFunctionReturn(0); 1645 } 1646 1647 static PetscErrorCode TSTrajectorySetFromOptions_Memory(PetscOptionItems *PetscOptionsObject,TSTrajectory tj) 1648 { 1649 TJScheduler *tjsch = (TJScheduler*)tj->data; 1650 PetscErrorCode ierr; 1651 1652 PetscFunctionBegin; 1653 ierr = PetscOptionsHead(PetscOptionsObject,"Memory based TS trajectory options");CHKERRQ(ierr); 1654 { 1655 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); 1656 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); 1657 ierr = PetscOptionsInt("-ts_trajectory_stride","Stride to save checkpoints to file","TSTrajectorySetStride_Memory",tjsch->stride,&tjsch->stride,NULL);CHKERRQ(ierr); 1658 #if defined(PETSC_HAVE_REVOLVE) 1659 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); 1660 #endif 1661 ierr = PetscOptionsBool("-ts_trajectory_save_stack","Save all stack to disk","TSTrajectorySetSaveStack",tjsch->save_stack,&tjsch->save_stack,NULL);CHKERRQ(ierr); 1662 ierr = PetscOptionsBool("-ts_trajectory_solution_only","Checkpoint solution only","TSTrajectorySetSolutionOnly",tjsch->stack.solution_only,&tjsch->stack.solution_only,NULL);CHKERRQ(ierr); 1663 ierr = PetscOptionsBool("-ts_trajectory_use_dram","Use DRAM for checkpointing","TSTrajectorySetUseDRAM",tjsch->stack.use_dram,&tjsch->stack.use_dram,NULL);CHKERRQ(ierr); 1664 } 1665 ierr = PetscOptionsTail();CHKERRQ(ierr); 1666 PetscFunctionReturn(0); 1667 } 1668 1669 static PetscErrorCode TSTrajectorySetUp_Memory(TSTrajectory tj,TS ts) 1670 { 1671 TJScheduler *tjsch = (TJScheduler*)tj->data; 1672 Stack *stack = &tjsch->stack; 1673 #if defined(PETSC_HAVE_REVOLVE) 1674 RevolveCTX *rctx,*rctx2; 1675 DiskStack *diskstack = &tjsch->diskstack; 1676 PetscInt diskblocks; 1677 #endif 1678 PetscInt numY; 1679 PetscBool flg; 1680 PetscErrorCode ierr; 1681 1682 PetscFunctionBegin; 1683 PetscStrcmp(((PetscObject)ts->adapt)->type_name,TSADAPTNONE,&flg); 1684 if (flg) tjsch->total_steps = PetscMin(ts->max_steps,(PetscInt)(PetscCeilReal((ts->max_time-ts->ptime)/ts->time_step))); /* fixed time step */ 1685 if (tjsch->max_cps_ram > 0) stack->stacksize = tjsch->max_cps_ram; 1686 1687 if (tjsch->stride > 1) { /* two level mode */ 1688 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."); 1689 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 */ 1690 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 */ 1691 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 */ 1692 } else { /* single level mode */ 1693 if (flg) { /* fixed time step */ 1694 if (tjsch->max_cps_ram >= tjsch->total_steps-1 || tjsch->max_cps_ram < 1) tjsch->stype = NONE; /* checkpoint all */ 1695 else tjsch->stype = (tjsch->max_cps_disk>1) ? REVOLVE_MULTISTAGE : REVOLVE_OFFLINE; 1696 } else tjsch->stype = NONE; /* checkpoint all for adaptive time step */ 1697 #if defined(PETSC_HAVE_REVOLVE) 1698 if (tjsch->use_online) tjsch->stype = REVOLVE_ONLINE; /* trick into online (for testing purpose only) */ 1699 #endif 1700 } 1701 1702 if (tjsch->stype > TWO_LEVEL_NOREVOLVE) { 1703 #ifndef PETSC_HAVE_REVOLVE 1704 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."); 1705 #else 1706 switch (tjsch->stype) { 1707 case TWO_LEVEL_REVOLVE: 1708 revolve_create_offline(tjsch->stride,tjsch->max_cps_ram); 1709 break; 1710 case TWO_LEVEL_TWO_REVOLVE: 1711 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. */ 1712 diskstack->stacksize = diskblocks; 1713 revolve_create_offline(tjsch->stride,tjsch->max_cps_ram); 1714 revolve2_create_offline((tjsch->total_steps+tjsch->stride-1)/tjsch->stride,diskblocks); 1715 ierr = PetscCalloc1(1,&rctx2);CHKERRQ(ierr); 1716 rctx2->snaps_in = diskblocks; 1717 rctx2->reverseonestep = PETSC_FALSE; 1718 rctx2->check = 0; 1719 rctx2->oldcapo = 0; 1720 rctx2->capo = 0; 1721 rctx2->info = 2; 1722 rctx2->fine = (tjsch->total_steps+tjsch->stride-1)/tjsch->stride; 1723 tjsch->rctx2 = rctx2; 1724 diskstack->top = -1; 1725 ierr = PetscMalloc1(diskstack->stacksize*sizeof(PetscInt),&diskstack->container);CHKERRQ(ierr); 1726 break; 1727 case REVOLVE_OFFLINE: 1728 revolve_create_offline(tjsch->total_steps,tjsch->max_cps_ram); 1729 break; 1730 case REVOLVE_ONLINE: 1731 stack->stacksize = tjsch->max_cps_ram; 1732 revolve_create_online(tjsch->max_cps_ram); 1733 break; 1734 case REVOLVE_MULTISTAGE: 1735 revolve_create_multistage(tjsch->total_steps,tjsch->max_cps_ram+tjsch->max_cps_disk,tjsch->max_cps_ram); 1736 break; 1737 default: 1738 break; 1739 } 1740 ierr = PetscCalloc1(1,&rctx);CHKERRQ(ierr); 1741 rctx->snaps_in = tjsch->max_cps_ram; /* for theta methods snaps_in=2*max_cps_ram */ 1742 rctx->reverseonestep = PETSC_FALSE; 1743 rctx->check = 0; 1744 rctx->oldcapo = 0; 1745 rctx->capo = 0; 1746 rctx->info = 2; 1747 rctx->fine = (tjsch->stride > 1) ? tjsch->stride : tjsch->total_steps; 1748 tjsch->rctx = rctx; 1749 if (tjsch->stype == REVOLVE_ONLINE) rctx->fine = -1; 1750 #endif 1751 } else { 1752 if (tjsch->stype == TWO_LEVEL_NOREVOLVE) stack->stacksize = tjsch->stride-1; /* need tjsch->stride-1 at most */ 1753 if (tjsch->stype == NONE) { 1754 if (flg) stack->stacksize = stack->solution_only ? tjsch->total_steps : tjsch->total_steps-1; /* fix time step */ 1755 else { /* adaptive time step */ 1756 /* if max_cps_ram is not specified, use maximal allowed number of steps for stack size */ 1757 if(tjsch->max_cps_ram == -1) stack->stacksize = ts->max_steps < PETSC_MAX_INT ? ts->max_steps : 10000; 1758 tjsch->total_steps = stack->solution_only ? stack->stacksize : stack->stacksize+1; /* will be updated as time integration advances */ 1759 } 1760 } 1761 } 1762 1763 tjsch->recompute = PETSC_FALSE; 1764 ierr = TSGetStages(ts,&numY,PETSC_IGNORE);CHKERRQ(ierr); 1765 ierr = StackCreate(stack,stack->stacksize,numY);CHKERRQ(ierr); 1766 PetscFunctionReturn(0); 1767 } 1768 1769 static PetscErrorCode TSTrajectoryDestroy_Memory(TSTrajectory tj) 1770 { 1771 TJScheduler *tjsch = (TJScheduler*)tj->data; 1772 PetscErrorCode ierr; 1773 1774 PetscFunctionBegin; 1775 if (tjsch->stype > TWO_LEVEL_NOREVOLVE) { 1776 #if defined(PETSC_HAVE_REVOLVE) 1777 revolve_reset(); 1778 if (tjsch->stype == TWO_LEVEL_TWO_REVOLVE) { 1779 revolve2_reset(); 1780 ierr = PetscFree(tjsch->diskstack.container);CHKERRQ(ierr); 1781 } 1782 #endif 1783 } 1784 ierr = StackDestroy(&tjsch->stack);CHKERRQ(ierr); 1785 #if defined(PETSC_HAVE_REVOLVE) 1786 if (tjsch->stype > TWO_LEVEL_NOREVOLVE) { 1787 ierr = PetscFree(tjsch->rctx);CHKERRQ(ierr); 1788 ierr = PetscFree(tjsch->rctx2);CHKERRQ(ierr); 1789 } 1790 #endif 1791 ierr = PetscFree(tjsch);CHKERRQ(ierr); 1792 PetscFunctionReturn(0); 1793 } 1794 1795 /*MC 1796 TSTRAJECTORYMEMORY - Stores each solution of the ODE/ADE in memory 1797 1798 Level: intermediate 1799 1800 .seealso: TSTrajectoryCreate(), TS, TSTrajectorySetType() 1801 1802 M*/ 1803 PETSC_EXTERN PetscErrorCode TSTrajectoryCreate_Memory(TSTrajectory tj,TS ts) 1804 { 1805 TJScheduler *tjsch; 1806 PetscErrorCode ierr; 1807 1808 PetscFunctionBegin; 1809 tj->ops->set = TSTrajectorySet_Memory; 1810 tj->ops->get = TSTrajectoryGet_Memory; 1811 tj->ops->setup = TSTrajectorySetUp_Memory; 1812 tj->ops->destroy = TSTrajectoryDestroy_Memory; 1813 tj->ops->setfromoptions = TSTrajectorySetFromOptions_Memory; 1814 1815 ierr = PetscCalloc1(1,&tjsch);CHKERRQ(ierr); 1816 tjsch->stype = NONE; 1817 tjsch->max_cps_ram = -1; /* -1 indicates that it is not set */ 1818 tjsch->max_cps_disk = -1; /* -1 indicates that it is not set */ 1819 tjsch->stride = 0; /* if not zero, two-level checkpointing will be used */ 1820 #if defined(PETSC_HAVE_REVOLVE) 1821 tjsch->use_online = PETSC_FALSE; 1822 #endif 1823 tjsch->save_stack = PETSC_TRUE; 1824 1825 tjsch->stack.solution_only = PETSC_TRUE; 1826 1827 tj->data = tjsch; 1828 PetscFunctionReturn(0); 1829 } 1830