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