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