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