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; 17 } *StackElement; 18 19 #ifdef PETSC_HAVE_REVOLVE 20 typedef struct _RevolveCTX { 21 PetscBool reverseonestep; 22 PetscInt where; 23 PetscInt snaps_in; 24 PetscInt stepsleft; 25 PetscInt check; 26 PetscInt oldcapo; 27 PetscInt capo; 28 PetscInt fine; 29 PetscInt info; 30 } RevolveCTX; 31 #endif 32 33 typedef struct _Stack { 34 PetscInt stacksize; 35 PetscInt top; 36 StackElement *container; 37 PetscInt numY; 38 PetscBool solution_only; 39 } Stack; 40 41 typedef struct _DiskStack { 42 PetscInt stacksize; 43 PetscInt top; 44 PetscInt *container; 45 } DiskStack; 46 47 typedef struct _TJScheduler { 48 SchedulerType stype; 49 #ifdef PETSC_HAVE_REVOLVE 50 RevolveCTX *rctx,*rctx2; 51 PetscBool use_online; 52 #endif 53 PetscBool recompute; 54 PetscBool skip_trajectory; 55 PetscBool save_stack; 56 MPI_Comm comm; 57 PetscInt max_cps_ram; /* maximum checkpoints in RAM */ 58 PetscInt max_cps_disk; /* maximum checkpoints on disk */ 59 PetscInt stride; 60 PetscInt total_steps; /* total number of steps */ 61 Stack stack; 62 DiskStack diskstack; 63 } TJScheduler; 64 65 #undef __FUNCT__ 66 #define __FUNCT__ "StackCreate" 67 static PetscErrorCode StackCreate(Stack *stack,PetscInt size,PetscInt ny) 68 { 69 PetscErrorCode ierr; 70 71 PetscFunctionBegin; 72 stack->top = -1; 73 stack->numY = ny; 74 75 ierr = PetscMalloc1(size*sizeof(StackElement),&stack->container);CHKERRQ(ierr); 76 PetscFunctionReturn(0); 77 } 78 79 #undef __FUNCT__ 80 #define __FUNCT__ "StackDestroy" 81 static PetscErrorCode StackDestroy(Stack *stack) 82 { 83 PetscInt i; 84 PetscErrorCode ierr; 85 86 PetscFunctionBegin; 87 if (stack->top > -1) { 88 for (i=0;i<=stack->top;i++) { 89 ierr = VecDestroy(&stack->container[i]->X);CHKERRQ(ierr); 90 if (!stack->solution_only) { 91 ierr = VecDestroyVecs(stack->numY,&stack->container[i]->Y);CHKERRQ(ierr); 92 } 93 ierr = PetscFree(stack->container[i]);CHKERRQ(ierr); 94 } 95 } 96 ierr = PetscFree(stack->container);CHKERRQ(ierr); 97 PetscFunctionReturn(0); 98 } 99 100 #undef __FUNCT__ 101 #define __FUNCT__ "StackPush" 102 static PetscErrorCode StackPush(Stack *stack,StackElement e) 103 { 104 PetscFunctionBegin; 105 if (stack->top+1 >= stack->stacksize) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_MEMC,"Maximum stack size (%D) exceeded",stack->stacksize); 106 stack->container[++stack->top] = e; 107 PetscFunctionReturn(0); 108 } 109 110 #undef __FUNCT__ 111 #define __FUNCT__ "StackPop" 112 static PetscErrorCode StackPop(Stack *stack,StackElement *e) 113 { 114 PetscFunctionBegin; 115 if (stack->top == -1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_MEMC,"Empty stack"); 116 *e = stack->container[stack->top--]; 117 PetscFunctionReturn(0); 118 } 119 120 #undef __FUNCT__ 121 #define __FUNCT__ "StackTop" 122 static PetscErrorCode StackTop(Stack *stack,StackElement *e) 123 { 124 PetscFunctionBegin; 125 *e = stack->container[stack->top]; 126 PetscFunctionReturn(0); 127 } 128 129 #ifdef PETSC_HAVE_REVOLVE 130 #undef __FUNCT__ 131 #define __FUNCT__ "StackFind" 132 static PetscErrorCode StackFind(Stack *stack,StackElement *e,PetscInt index) 133 { 134 PetscFunctionBegin; 135 *e = stack->container[index]; 136 PetscFunctionReturn(0); 137 } 138 #endif 139 140 #undef __FUNCT__ 141 #define __FUNCT__ "OutputBIN" 142 static PetscErrorCode OutputBIN(const char *filename,PetscViewer *viewer) 143 { 144 PetscErrorCode ierr; 145 146 PetscFunctionBegin; 147 ierr = PetscViewerCreate(PETSC_COMM_WORLD,viewer);CHKERRQ(ierr); 148 ierr = PetscViewerSetType(*viewer,PETSCVIEWERBINARY);CHKERRQ(ierr); 149 ierr = PetscViewerFileSetMode(*viewer,FILE_MODE_WRITE);CHKERRQ(ierr); 150 ierr = PetscViewerFileSetName(*viewer,filename);CHKERRQ(ierr); 151 PetscFunctionReturn(0); 152 } 153 154 #undef __FUNCT__ 155 #define __FUNCT__ "WriteToDisk" 156 static PetscErrorCode WriteToDisk(PetscInt stepnum,PetscReal time,PetscReal timeprev,Vec X,Vec *Y,PetscInt numY,PetscBool solution_only,PetscViewer viewer) 157 { 158 PetscInt i; 159 PetscErrorCode ierr; 160 161 PetscFunctionBegin; 162 ierr = PetscViewerBinaryWrite(viewer,&stepnum,1,PETSC_INT,PETSC_FALSE);CHKERRQ(ierr); 163 ierr = VecView(X,viewer);CHKERRQ(ierr); 164 for (i=0;!solution_only && i<numY;i++) { 165 ierr = VecView(Y[i],viewer);CHKERRQ(ierr); 166 } 167 ierr = PetscViewerBinaryWrite(viewer,&time,1,PETSC_REAL,PETSC_FALSE);CHKERRQ(ierr); 168 ierr = PetscViewerBinaryWrite(viewer,&timeprev,1,PETSC_REAL,PETSC_FALSE);CHKERRQ(ierr); 169 PetscFunctionReturn(0); 170 } 171 172 #undef __FUNCT__ 173 #define __FUNCT__ "ReadFromDisk" 174 static PetscErrorCode ReadFromDisk(PetscInt *stepnum,PetscReal *time,PetscReal *timeprev,Vec X,Vec *Y,PetscInt numY,PetscBool solution_only,PetscViewer viewer) 175 { 176 PetscInt i; 177 PetscErrorCode ierr; 178 179 PetscFunctionBegin; 180 ierr = PetscViewerBinaryRead(viewer,stepnum,1,NULL,PETSC_INT);CHKERRQ(ierr); 181 ierr = VecLoad(X,viewer);CHKERRQ(ierr); 182 for (i=0;!solution_only && i<numY;i++) { 183 ierr = VecLoad(Y[i],viewer);CHKERRQ(ierr); 184 } 185 ierr = PetscViewerBinaryRead(viewer,time,1,NULL,PETSC_REAL);CHKERRQ(ierr); 186 ierr = PetscViewerBinaryRead(viewer,timeprev,1,NULL,PETSC_REAL);CHKERRQ(ierr); 187 PetscFunctionReturn(0); 188 } 189 190 #undef __FUNCT__ 191 #define __FUNCT__ "StackDumpAll" 192 static PetscErrorCode StackDumpAll(TS ts,Stack *stack,PetscInt id) 193 { 194 Vec *Y; 195 PetscInt i; 196 StackElement e; 197 PetscViewer viewer; 198 char filename[PETSC_MAX_PATH_LEN]; 199 PetscErrorCode ierr; 200 201 PetscFunctionBegin; 202 if (id == 1) { 203 PetscMPIInt rank; 204 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)ts),&rank);CHKERRQ(ierr); 205 if (!rank) { 206 ierr = PetscRMTree("SA-data");CHKERRQ(ierr); 207 ierr = PetscMkdir("SA-data");CHKERRQ(ierr); 208 } 209 } 210 ierr = PetscSNPrintf(filename,sizeof(filename),"SA-data/SA-STACK%06d.bin",id);CHKERRQ(ierr); 211 ierr = OutputBIN(filename,&viewer);CHKERRQ(ierr); 212 for (i=0;i<stack->stacksize;i++) { 213 e = stack->container[i]; 214 ierr = WriteToDisk(e->stepnum,e->time,e->timeprev,e->X,e->Y,stack->numY,stack->solution_only,viewer);CHKERRQ(ierr); 215 } 216 /* 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 */ 217 ierr = TSGetStages(ts,&stack->numY,&Y);CHKERRQ(ierr); 218 ierr = WriteToDisk(ts->total_steps,ts->ptime,ts->ptime_prev,ts->vec_sol,Y,stack->numY,stack->solution_only,viewer);CHKERRQ(ierr); 219 for (i=0;i<stack->stacksize;i++) { 220 ierr = StackPop(stack,&e);CHKERRQ(ierr); 221 ierr = VecDestroy(&e->X);CHKERRQ(ierr); 222 if (!stack->solution_only) { 223 ierr = VecDestroyVecs(stack->numY,&e->Y);CHKERRQ(ierr); 224 } 225 ierr = PetscFree(e);CHKERRQ(ierr); 226 } 227 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 228 PetscFunctionReturn(0); 229 } 230 231 #undef __FUNCT__ 232 #define __FUNCT__ "StackLoadAll" 233 static PetscErrorCode StackLoadAll(TS ts,Stack *stack,PetscInt id) 234 { 235 Vec *Y; 236 PetscInt i; 237 StackElement e; 238 PetscViewer viewer; 239 char filename[PETSC_MAX_PATH_LEN]; 240 PetscErrorCode ierr; 241 242 PetscFunctionBegin; 243 ierr = PetscSNPrintf(filename,sizeof filename,"SA-data/SA-STACK%06d.bin",id);CHKERRQ(ierr); 244 ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);CHKERRQ(ierr); 245 for (i=0;i<stack->stacksize;i++) { 246 ierr = PetscCalloc1(1,&e); 247 ierr = TSGetStages(ts,&stack->numY,&Y);CHKERRQ(ierr); 248 ierr = VecDuplicate(Y[0],&e->X);CHKERRQ(ierr); 249 if (!stack->solution_only && stack->numY>0) { 250 ierr = VecDuplicateVecs(Y[0],stack->numY,&e->Y);CHKERRQ(ierr); 251 } 252 ierr = StackPush(stack,e);CHKERRQ(ierr); 253 } 254 for (i=0;i<stack->stacksize;i++) { 255 e = stack->container[i]; 256 ierr = ReadFromDisk(&e->stepnum,&e->time,&e->timeprev,e->X,e->Y,stack->numY,stack->solution_only,viewer);CHKERRQ(ierr); 257 } 258 /* load the last step into TS */ 259 ierr = TSGetStages(ts,&stack->numY,&Y);CHKERRQ(ierr); 260 ierr = ReadFromDisk(&ts->total_steps,&ts->ptime,&ts->ptime_prev,ts->vec_sol,Y,stack->numY,stack->solution_only,viewer);CHKERRQ(ierr); 261 ierr = TSSetTimeStep(ts,ts->ptime_prev-ts->ptime);CHKERRQ(ierr); 262 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 263 PetscFunctionReturn(0); 264 } 265 266 #undef __FUNCT__ 267 #define __FUNCT__ "DumpSingle" 268 static PetscErrorCode DumpSingle(TS ts,Stack *stack,PetscInt id) 269 { 270 Vec *Y; 271 PetscInt stepnum; 272 PetscViewer viewer; 273 char filename[PETSC_MAX_PATH_LEN]; 274 PetscErrorCode ierr; 275 276 PetscFunctionBegin; 277 ierr = TSGetTotalSteps(ts,&stepnum);CHKERRQ(ierr); 278 if (id == 1) { 279 PetscMPIInt rank; 280 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)ts),&rank);CHKERRQ(ierr); 281 if (!rank) { 282 ierr = PetscRMTree("SA-data");CHKERRQ(ierr); 283 ierr = PetscMkdir("SA-data");CHKERRQ(ierr); 284 } 285 } 286 ierr = PetscSNPrintf(filename,sizeof(filename),"SA-data/SA-CPS%06d.bin",id);CHKERRQ(ierr); 287 ierr = OutputBIN(filename,&viewer);CHKERRQ(ierr); 288 289 ierr = PetscLogEventBegin(Disk_Write,ts,0,0,0);CHKERRQ(ierr); 290 ierr = TSGetStages(ts,&stack->numY,&Y);CHKERRQ(ierr); 291 ierr = WriteToDisk(stepnum,ts->ptime,ts->ptime_prev,ts->vec_sol,Y,stack->numY,stack->solution_only,viewer);CHKERRQ(ierr); 292 ierr = PetscLogEventEnd(Disk_Write,ts,0,0,0);CHKERRQ(ierr); 293 294 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 295 PetscFunctionReturn(0); 296 } 297 298 #undef __FUNCT__ 299 #define __FUNCT__ "LoadSingle" 300 static PetscErrorCode LoadSingle(TS ts,Stack *stack,PetscInt id) 301 { 302 Vec *Y; 303 PetscViewer viewer; 304 char filename[PETSC_MAX_PATH_LEN]; 305 PetscErrorCode ierr; 306 307 PetscFunctionBegin; 308 ierr = PetscSNPrintf(filename,sizeof filename,"SA-data/SA-CPS%06d.bin",id);CHKERRQ(ierr); 309 ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);CHKERRQ(ierr); 310 311 ierr = PetscLogEventBegin(Disk_Read,ts,0,0,0);CHKERRQ(ierr); 312 ierr = TSGetStages(ts,&stack->numY,&Y);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(Disk_Read,ts,0,0,0);CHKERRQ(ierr); 315 316 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 317 PetscFunctionReturn(0); 318 } 319 320 #undef __FUNCT__ 321 #define __FUNCT__ "ElementCreate" 322 static PetscErrorCode ElementCreate(TS ts,Stack *stack,StackElement *e,PetscInt stepnum,PetscReal time,Vec X) 323 { 324 Vec *Y; 325 PetscInt i; 326 PetscReal timeprev; 327 PetscErrorCode ierr; 328 329 PetscFunctionBegin; 330 ierr = PetscCalloc1(1,e);CHKERRQ(ierr); 331 ierr = VecDuplicate(X,&(*e)->X);CHKERRQ(ierr); 332 ierr = VecCopy(X,(*e)->X);CHKERRQ(ierr); 333 if (stack->numY > 0 && !stack->solution_only) { 334 ierr = TSGetStages(ts,&stack->numY,&Y);CHKERRQ(ierr); 335 ierr = VecDuplicateVecs(Y[0],stack->numY,&(*e)->Y);CHKERRQ(ierr); 336 for (i=0;i<stack->numY;i++) { 337 ierr = VecCopy(Y[i],(*e)->Y[i]);CHKERRQ(ierr); 338 } 339 } 340 (*e)->stepnum = stepnum; 341 (*e)->time = time; 342 /* for consistency */ 343 if (stepnum == 0) { 344 (*e)->timeprev = (*e)->time - ts->time_step; 345 } else { 346 ierr = TSGetPrevTime(ts,&timeprev);CHKERRQ(ierr); 347 (*e)->timeprev = timeprev; 348 } 349 PetscFunctionReturn(0); 350 } 351 352 #undef __FUNCT__ 353 #define __FUNCT__ "ElementDestroy" 354 static PetscErrorCode ElementDestroy(Stack *stack,StackElement e) 355 { 356 PetscErrorCode ierr; 357 358 PetscFunctionBegin; 359 ierr = VecDestroy(&e->X);CHKERRQ(ierr); 360 if (!stack->solution_only) { 361 ierr = VecDestroyVecs(stack->numY,&e->Y);CHKERRQ(ierr); 362 } 363 ierr = PetscFree(e);CHKERRQ(ierr); 364 PetscFunctionReturn(0); 365 } 366 367 #undef __FUNCT__ 368 #define __FUNCT__ "UpdateTS" 369 static PetscErrorCode UpdateTS(TS ts,Stack *stack,StackElement e) 370 { 371 Vec *Y; 372 PetscInt i; 373 PetscErrorCode ierr; 374 375 PetscFunctionBegin; 376 ierr = VecCopy(e->X,ts->vec_sol);CHKERRQ(ierr); 377 if (!stack->solution_only) { 378 ierr = TSGetStages(ts,&stack->numY,&Y);CHKERRQ(ierr); 379 for (i=0;i<stack->numY;i++) { 380 ierr = VecCopy(e->Y[i],Y[i]);CHKERRQ(ierr); 381 } 382 } 383 ierr = TSSetTimeStep(ts,e->timeprev-e->time);CHKERRQ(ierr); /* stepsize will be negative */ 384 ts->ptime = e->time; 385 ts->ptime_prev = e->timeprev; 386 PetscFunctionReturn(0); 387 } 388 389 #undef __FUNCT__ 390 #define __FUNCT__ "ReCompute" 391 static PetscErrorCode ReCompute(TS ts,TJScheduler *tjsch,PetscInt stepnumbegin,PetscInt stepnumend) 392 { 393 Stack *stack = &tjsch->stack; 394 PetscInt i,adjsteps; 395 PetscReal stepsize; 396 PetscErrorCode ierr; 397 398 PetscFunctionBegin; 399 adjsteps = ts->steps; 400 /* reset ts context */ 401 ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr); 402 ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr); 403 ts->steps = stepnumbegin; /* global step number */ 404 for (i=ts->steps;i<stepnumend;i++) { /* assume fixed step size */ 405 if (stack->solution_only && !tjsch->skip_trajectory) { /* revolve online need this */ 406 ierr = TSTrajectorySet(ts->trajectory,ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 407 } 408 ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 409 ierr = TSStep(ts);CHKERRQ(ierr); 410 if (!stack->solution_only && !tjsch->skip_trajectory) { 411 ierr = TSTrajectorySet(ts->trajectory,ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 412 } 413 if (ts->event) { 414 ierr = TSEventMonitor(ts);CHKERRQ(ierr); 415 } 416 if (!ts->steprollback) { 417 ierr = TSPostStep(ts);CHKERRQ(ierr); 418 } 419 } 420 ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr); 421 ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr); 422 ts->steps = adjsteps; 423 ts->total_steps = stepnumend; 424 PetscFunctionReturn(0); 425 } 426 427 #undef __FUNCT__ 428 #define __FUNCT__ "TopLevelStoreAll" 429 static PetscErrorCode TopLevelStoreAll(TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscInt localstepnum,PetscInt laststridesize,PetscBool *reset) 430 { 431 Stack *stack = &tjsch->stack; 432 PetscInt stridenum; 433 PetscErrorCode ierr; 434 435 PetscFunctionBegin; 436 stridenum = stepnum/tjsch->stride; 437 /** 438 * make sure saved checkpoint id starts from 1 439 * skip last stride when using stridenum+1 440 * skip first stride when using stridenum 441 **/ 442 if (stack->solution_only) { 443 if (tjsch->save_stack) { 444 if (localstepnum == tjsch->stride-1 && stepnum < tjsch->total_steps-laststridesize) { /* current step will be saved without going through stack */ 445 ierr = StackDumpAll(ts,stack,stridenum+1);CHKERRQ(ierr); 446 if (reset) *reset = PETSC_TRUE; 447 } 448 } else { 449 if (!tjsch->recompute && localstepnum == 0 && stepnum < tjsch->total_steps-laststridesize) { 450 ierr = DumpSingle(ts,stack,stridenum+1);CHKERRQ(ierr); 451 } 452 } 453 } else { 454 if (tjsch->save_stack) { 455 if (localstepnum == 0 && stepnum < tjsch->total_steps) { /* first stride skipped outside of this function */ 456 ierr = StackDumpAll(ts,stack,stridenum);CHKERRQ(ierr); 457 if (reset) *reset = PETSC_TRUE; 458 } 459 } else { 460 if (!tjsch->recompute && localstepnum == 1 && stepnum < tjsch->total_steps-laststridesize) { 461 ierr = DumpSingle(ts,stack,stridenum+1);CHKERRQ(ierr); 462 } 463 } 464 } 465 PetscFunctionReturn(0); 466 } 467 468 #undef __FUNCT__ 469 #define __FUNCT__ "SetTrajN" 470 static PetscErrorCode SetTrajN(TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X) 471 { 472 Stack *stack = &tjsch->stack; 473 StackElement e; 474 PetscErrorCode ierr; 475 476 PetscFunctionBegin; 477 /* skip the last two steps of each stride or the whole interval */ 478 if (stack->solution_only && (stepnum >= tjsch->total_steps-1 || tjsch->recompute)) PetscFunctionReturn(0); //? 479 /* skip the first and the last steps of each stride or the whole interval */ 480 if (!stack->solution_only && (stepnum == 0 || stepnum == tjsch->total_steps)) PetscFunctionReturn(0); 481 482 if (stepnum < stack->top) SETERRQ(tjsch->comm,PETSC_ERR_MEMC,"Illegal modification of a non-top stack element"); 483 ierr = ElementCreate(ts,stack,&e,stepnum,time,X);CHKERRQ(ierr); 484 ierr = StackPush(stack,e);CHKERRQ(ierr); 485 PetscFunctionReturn(0); 486 } 487 488 #undef __FUNCT__ 489 #define __FUNCT__ "GetTrajN" 490 static PetscErrorCode GetTrajN(TS ts,TJScheduler *tjsch,PetscInt stepnum) 491 { 492 Stack *stack = &tjsch->stack; 493 PetscReal stepsize; 494 StackElement e; 495 PetscErrorCode ierr; 496 497 PetscFunctionBegin; 498 if (stepnum == tjsch->total_steps) { 499 ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr); 500 ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr); 501 PetscFunctionReturn(0); 502 } 503 /* restore a checkpoint */ 504 ierr = StackTop(stack,&e);CHKERRQ(ierr); 505 ierr = UpdateTS(ts,stack,e);CHKERRQ(ierr); 506 if (stack->solution_only) {/* recompute one step */ 507 tjsch->recompute = PETSC_TRUE; 508 ierr = ReCompute(ts,tjsch,e->stepnum,stepnum);CHKERRQ(ierr); 509 } 510 ierr = StackPop(stack,&e);CHKERRQ(ierr); 511 ierr = ElementDestroy(stack,e);CHKERRQ(ierr); 512 PetscFunctionReturn(0); 513 } 514 515 #undef __FUNCT__ 516 #define __FUNCT__ "SetTrajTLNR" 517 static PetscErrorCode SetTrajTLNR(TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X) 518 { 519 Stack *stack = &tjsch->stack; 520 PetscInt localstepnum,laststridesize; 521 StackElement e; 522 PetscErrorCode ierr; 523 524 PetscFunctionBegin; 525 if (!stack->solution_only && stepnum == 0) PetscFunctionReturn(0); 526 if (stack->solution_only && stepnum == tjsch->total_steps) PetscFunctionReturn(0); 527 if (tjsch->save_stack && tjsch->recompute) PetscFunctionReturn(0); 528 529 localstepnum = stepnum%tjsch->stride; 530 /* (stride size-1) checkpoints are saved in each stride; an extra point is added by StackDumpAll() */ 531 laststridesize = tjsch->total_steps%tjsch->stride; 532 if (!laststridesize) laststridesize = tjsch->stride; 533 534 if (!tjsch->recompute) { 535 TopLevelStoreAll(ts,tjsch,stepnum,localstepnum,laststridesize,PETSC_NULL); 536 if (!tjsch->save_stack && stepnum < tjsch->total_steps-laststridesize) PetscFunctionReturn(0); 537 } 538 if (!stack->solution_only && localstepnum == 0) PetscFunctionReturn(0); /* skip last point in each stride at recompute stage or last stride */ 539 if (stack->solution_only && localstepnum == tjsch->stride-1) PetscFunctionReturn(0); /* skip last step in each stride at recompute stage or last stride */ 540 541 ierr = ElementCreate(ts,stack,&e,stepnum,time,X);CHKERRQ(ierr); 542 ierr = StackPush(stack,e);CHKERRQ(ierr); 543 PetscFunctionReturn(0); 544 } 545 546 #undef __FUNCT__ 547 #define __FUNCT__ "GetTrajTLNR" 548 static PetscErrorCode GetTrajTLNR(TS ts,TJScheduler *tjsch,PetscInt stepnum) 549 { 550 Stack *stack = &tjsch->stack; 551 PetscInt id,localstepnum,laststridesize; 552 PetscReal stepsize; 553 StackElement e; 554 PetscErrorCode ierr; 555 556 PetscFunctionBegin; 557 if (stepnum == tjsch->total_steps) { 558 ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr); 559 ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr); 560 PetscFunctionReturn(0); 561 } 562 563 localstepnum = stepnum%tjsch->stride; 564 laststridesize = tjsch->total_steps%tjsch->stride; 565 if (!laststridesize) laststridesize = tjsch->stride; 566 if (stack->solution_only) { 567 /* fill stack with info */ 568 if (localstepnum == 0 && tjsch->total_steps-stepnum >= laststridesize) { 569 id = stepnum/tjsch->stride; 570 if (tjsch->save_stack) { 571 ierr = StackLoadAll(ts,stack,id);CHKERRQ(ierr); 572 tjsch->recompute = PETSC_TRUE; 573 tjsch->skip_trajectory = PETSC_TRUE; 574 ierr = ReCompute(ts,tjsch,id*tjsch->stride-1,id*tjsch->stride);CHKERRQ(ierr); 575 tjsch->skip_trajectory = PETSC_FALSE; 576 } else { 577 ierr = LoadSingle(ts,stack,id);CHKERRQ(ierr); 578 tjsch->recompute = PETSC_TRUE; 579 ierr = ReCompute(ts,tjsch,(id-1)*tjsch->stride,id*tjsch->stride);CHKERRQ(ierr); 580 } 581 PetscFunctionReturn(0); 582 } 583 /* restore a checkpoint */ 584 ierr = StackPop(stack,&e);CHKERRQ(ierr); 585 ierr = UpdateTS(ts,stack,e);CHKERRQ(ierr); 586 tjsch->recompute = PETSC_TRUE; 587 tjsch->skip_trajectory = PETSC_TRUE; 588 ierr = ReCompute(ts,tjsch,e->stepnum,stepnum);CHKERRQ(ierr); 589 tjsch->skip_trajectory = PETSC_FALSE; 590 ierr = ElementDestroy(stack,e);CHKERRQ(ierr); 591 } else { 592 /* fill stack with info */ 593 if (localstepnum == 0 && tjsch->total_steps-stepnum >= laststridesize) { 594 id = stepnum/tjsch->stride; 595 if (tjsch->save_stack) { 596 ierr = StackLoadAll(ts,stack,id);CHKERRQ(ierr); 597 } else { 598 ierr = LoadSingle(ts,stack,id);CHKERRQ(ierr); 599 ierr = ElementCreate(ts,stack,&e,(id-1)*tjsch->stride+1,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 600 ierr = StackPush(stack,e);CHKERRQ(ierr); 601 tjsch->recompute = PETSC_TRUE; 602 ierr = ReCompute(ts,tjsch,e->stepnum,id*tjsch->stride);CHKERRQ(ierr); 603 } 604 PetscFunctionReturn(0); 605 } 606 /* restore a checkpoint */ 607 ierr = StackPop(stack,&e);CHKERRQ(ierr); 608 ierr = UpdateTS(ts,stack,e);CHKERRQ(ierr); 609 ierr = ElementDestroy(stack,e);CHKERRQ(ierr); 610 } 611 PetscFunctionReturn(0); 612 } 613 614 #ifdef PETSC_HAVE_REVOLVE 615 static void printwhattodo(PetscInt whattodo,RevolveCTX *rctx,PetscInt shift) 616 { 617 switch(whattodo) { 618 case 1: 619 PetscPrintf(PETSC_COMM_WORLD,"\x1B[35mAdvance from %D to %D\033[0m\n",rctx->oldcapo+shift,rctx->capo+shift); 620 break; 621 case 2: 622 PetscPrintf(PETSC_COMM_WORLD,"\x1B[35mStore in checkpoint number %D (located in RAM)\033[0m\n",rctx->check); 623 break; 624 case 3: 625 PetscPrintf(PETSC_COMM_WORLD,"\x1B[35mFirst turn: Initialize adjoints and reverse first step\033[0m\n"); 626 break; 627 case 4: 628 PetscPrintf(PETSC_COMM_WORLD,"\x1B[35mForward and reverse one step\033[0m\n"); 629 break; 630 case 5: 631 PetscPrintf(PETSC_COMM_WORLD,"\x1B[35mRestore in checkpoint number %D (located in RAM)\033[0m\n",rctx->check); 632 break; 633 case 7: 634 PetscPrintf(PETSC_COMM_WORLD,"\x1B[35mStore in checkpoint number %D (located on disk)\033[0m\n",rctx->check); 635 break; 636 case 8: 637 PetscPrintf(PETSC_COMM_WORLD,"\x1B[35mRestore in checkpoint number %D (located on disk)\033[0m\n",rctx->check); 638 break; 639 case -1: 640 PetscPrintf(PETSC_COMM_WORLD,"\x1B[35mError!"); 641 break; 642 } 643 } 644 645 #undef __FUNCT__ 646 #define __FUNCT__ "ApplyRevolve" 647 static PetscErrorCode ApplyRevolve(SchedulerType stype,RevolveCTX *rctx,PetscInt total_steps,PetscInt stepnum,PetscInt localstepnum,PetscInt *store) 648 { 649 PetscInt shift,whattodo; 650 651 PetscFunctionBegin; 652 *store = 0; 653 // if (rctx->reverseonestep && stepnum == total_steps) { /* intermediate information is ready inside TS, this happens at last time step */ 654 // rctx->reverseonestep = PETSC_FALSE; 655 // PetscFunctionReturn(0); 656 // } 657 if (rctx->stepsleft > 0) { /* advance the solution without checkpointing anything as Revolve requires */ 658 rctx->stepsleft--; 659 PetscFunctionReturn(0); 660 } 661 /* let Revolve determine what to do next */ 662 shift = stepnum-localstepnum; 663 rctx->oldcapo = rctx->capo; 664 rctx->capo = localstepnum; 665 whattodo = revolve_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where); 666 if (stype == REVOLVE_ONLINE && whattodo == 8) whattodo = 5; 667 if (stype == REVOLVE_ONLINE && whattodo == 7) whattodo = 2; 668 printwhattodo(whattodo,rctx,shift); 669 if (whattodo == -1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_LIB,"Error in the Revolve library"); 670 if (whattodo == 1) { /* advance some time steps */ 671 if (stype == REVOLVE_ONLINE && rctx->capo >= total_steps-1) { 672 revolve_turn(total_steps,&rctx->capo,&rctx->fine); 673 printwhattodo(whattodo,rctx,shift); 674 } 675 rctx->stepsleft = rctx->capo-rctx->oldcapo-1; 676 } 677 if (whattodo == 3 || whattodo == 4) { /* ready for a reverse step */ 678 rctx->reverseonestep = PETSC_TRUE; 679 } 680 if (whattodo == 5) { /* restore a checkpoint and ask Revolve what to do next */ 681 rctx->oldcapo = rctx->capo; 682 whattodo = revolve_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where); /* must return 1 or 3 or 4*/ 683 printwhattodo(whattodo,rctx,shift); 684 if (whattodo == 3 || whattodo == 4) rctx->reverseonestep = PETSC_TRUE; 685 if (whattodo == 1) rctx->stepsleft = rctx->capo-rctx->oldcapo; 686 } 687 if (whattodo == 7) { /* save the checkpoint to disk */ 688 *store = 2; 689 rctx->oldcapo = rctx->capo; 690 whattodo = revolve_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where); /* must return 1 */ 691 printwhattodo(whattodo,rctx,shift); 692 rctx->stepsleft = rctx->capo-rctx->oldcapo-1; 693 } 694 if (whattodo == 2) { /* store a checkpoint to RAM and ask Revolve how many time steps to advance next */ 695 *store = 1; 696 rctx->oldcapo = rctx->capo; 697 whattodo = revolve_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where); /* must return 1 */ 698 printwhattodo(whattodo,rctx,shift); 699 if (stype == REVOLVE_ONLINE && rctx->capo >= total_steps-1) { 700 revolve_turn(total_steps,&rctx->capo,&rctx->fine); 701 printwhattodo(whattodo,rctx,shift); 702 } 703 rctx->stepsleft = rctx->capo-rctx->oldcapo-1; 704 } 705 PetscFunctionReturn(0); 706 } 707 708 #undef __FUNCT__ 709 #define __FUNCT__ "SetTrajROF" 710 static PetscErrorCode SetTrajROF(TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X) 711 { 712 Stack *stack = &tjsch->stack; 713 PetscInt store; 714 StackElement e; 715 PetscErrorCode ierr; 716 717 PetscFunctionBegin; 718 if (!stack->solution_only && stepnum == 0) PetscFunctionReturn(0); 719 if (stack->solution_only && stepnum == tjsch->total_steps) PetscFunctionReturn(0); 720 ierr = ApplyRevolve(tjsch->stype,tjsch->rctx,tjsch->total_steps,stepnum,stepnum,&store);CHKERRQ(ierr); 721 if (store == 1) { 722 if (stepnum < stack->top) SETERRQ(tjsch->comm,PETSC_ERR_MEMC,"Illegal modification of a non-top stack element"); 723 ierr = ElementCreate(ts,stack,&e,stepnum,time,X);CHKERRQ(ierr); 724 ierr = StackPush(stack,e);CHKERRQ(ierr); 725 } 726 PetscFunctionReturn(0); 727 } 728 729 #undef __FUNCT__ 730 #define __FUNCT__ "GetTrajROF" 731 static PetscErrorCode GetTrajROF(TS ts,TJScheduler *tjsch,PetscInt stepnum) 732 { 733 Stack *stack = &tjsch->stack; 734 PetscInt whattodo,shift,store; 735 PetscReal stepsize; 736 StackElement e; 737 PetscErrorCode ierr; 738 739 PetscFunctionBegin; 740 if (stepnum == 0 || stepnum == tjsch->total_steps) { 741 ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr); 742 ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr); 743 if (tjsch->rctx->reverseonestep) tjsch->rctx->reverseonestep = PETSC_FALSE; 744 PetscFunctionReturn(0); 745 } 746 /* restore a checkpoint */ 747 ierr = StackTop(stack,&e);CHKERRQ(ierr); 748 ierr = UpdateTS(ts,stack,e);CHKERRQ(ierr); 749 if (stack->solution_only) { /* start with restoring a checkpoint */ 750 tjsch->rctx->capo = stepnum; 751 tjsch->rctx->oldcapo = tjsch->rctx->capo; 752 shift = 0; 753 whattodo = revolve_action(&tjsch->rctx->check,&tjsch->rctx->capo,&tjsch->rctx->fine,tjsch->rctx->snaps_in,&tjsch->rctx->info,&tjsch->rctx->where); 754 printwhattodo(whattodo,tjsch->rctx,shift); 755 } else { /* 2 revolve actions: restore a checkpoint and then advance */ 756 ierr = ApplyRevolve(tjsch->stype,tjsch->rctx,tjsch->total_steps,stepnum,stepnum,&store);CHKERRQ(ierr); 757 if (!tjsch->rctx->reverseonestep && tjsch->rctx->stepsleft > 0) { 758 tjsch->rctx->stepsleft--; 759 PetscPrintf(PETSC_COMM_WORLD,"\x1B[35mSkip the step from %D to %D (already checkpointed)\033[0m\n",tjsch->rctx->oldcapo,tjsch->rctx->oldcapo+1); 760 } 761 } 762 if (stack->solution_only || (!stack->solution_only && e->stepnum < stepnum)) { 763 tjsch->recompute = PETSC_TRUE; 764 ierr = ReCompute(ts,tjsch,e->stepnum,stepnum);CHKERRQ(ierr); 765 } 766 if ((stack->solution_only && e->stepnum+1 == stepnum) || (!stack->solution_only && e->stepnum == stepnum)) { 767 ierr = StackPop(stack,&e);CHKERRQ(ierr); 768 ierr = ElementDestroy(stack,e);CHKERRQ(ierr); 769 } 770 if (tjsch->rctx->reverseonestep) tjsch->rctx->reverseonestep = PETSC_FALSE; 771 PetscFunctionReturn(0); 772 } 773 774 #undef __FUNCT__ 775 #define __FUNCT__ "SetTrajRON" 776 static PetscErrorCode SetTrajRON(TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X) 777 { 778 Stack *stack = &tjsch->stack; 779 Vec *Y; 780 PetscInt i,store; 781 PetscReal timeprev; 782 StackElement e; 783 RevolveCTX *rctx = tjsch->rctx; 784 PetscErrorCode ierr; 785 786 PetscFunctionBegin; 787 if (!stack->solution_only && stepnum == 0) PetscFunctionReturn(0); 788 if (stack->solution_only && stepnum == tjsch->total_steps) PetscFunctionReturn(0); 789 ierr = ApplyRevolve(tjsch->stype,rctx,tjsch->total_steps,stepnum,stepnum,&store);CHKERRQ(ierr); 790 if (store == 1) { 791 if (rctx->check != stack->top+1) { /* overwrite some non-top checkpoint in the stack */ 792 ierr = StackFind(stack,&e,rctx->check);CHKERRQ(ierr); 793 ierr = VecCopy(X,e->X);CHKERRQ(ierr); 794 if (stack->numY > 0 && !stack->solution_only) { 795 ierr = TSGetStages(ts,&stack->numY,&Y);CHKERRQ(ierr); 796 for (i=0;i<stack->numY;i++) { 797 ierr = VecCopy(Y[i],e->Y[i]);CHKERRQ(ierr); 798 } 799 } 800 e->stepnum = stepnum; 801 e->time = time; 802 ierr = TSGetPrevTime(ts,&timeprev);CHKERRQ(ierr); 803 e->timeprev = timeprev; 804 } else { 805 if (stepnum < stack->top) SETERRQ(tjsch->comm,PETSC_ERR_MEMC,"Illegal modification of a non-top stack element"); 806 ierr = ElementCreate(ts,stack,&e,stepnum,time,X);CHKERRQ(ierr); 807 ierr = StackPush(stack,e);CHKERRQ(ierr); 808 } 809 } 810 PetscFunctionReturn(0); 811 } 812 813 #undef __FUNCT__ 814 #define __FUNCT__ "GetTrajRON" 815 static PetscErrorCode GetTrajRON(TS ts,TJScheduler *tjsch,PetscInt stepnum) 816 { 817 Stack *stack = &tjsch->stack; 818 PetscInt whattodo,shift; 819 PetscReal stepsize; 820 StackElement e; 821 PetscErrorCode ierr; 822 823 PetscFunctionBegin; 824 if (stepnum == 0 || stepnum == tjsch->total_steps) { 825 ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr); 826 ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr); 827 if (tjsch->rctx->reverseonestep) tjsch->rctx->reverseonestep = PETSC_FALSE; 828 PetscFunctionReturn(0); 829 } 830 tjsch->rctx->capo = stepnum; 831 tjsch->rctx->oldcapo = tjsch->rctx->capo; 832 shift = 0; 833 whattodo = revolve_action(&tjsch->rctx->check,&tjsch->rctx->capo,&tjsch->rctx->fine,tjsch->rctx->snaps_in,&tjsch->rctx->info,&tjsch->rctx->where); /* whattodo=restore */ 834 if (whattodo == 8) whattodo = 5; 835 printwhattodo(whattodo,tjsch->rctx,shift); 836 /* restore a checkpoint */ 837 ierr = StackFind(stack,&e,tjsch->rctx->check);CHKERRQ(ierr); 838 ierr = UpdateTS(ts,stack,e);CHKERRQ(ierr); 839 if (!stack->solution_only) { /* whattodo must be 5 */ 840 /* ask Revolve what to do next */ 841 tjsch->rctx->oldcapo = tjsch->rctx->capo; 842 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*/ 843 printwhattodo(whattodo,tjsch->rctx,shift); 844 if (whattodo == 3 || whattodo == 4) tjsch->rctx->reverseonestep = PETSC_TRUE; 845 if (whattodo == 1) tjsch->rctx->stepsleft = tjsch->rctx->capo-tjsch->rctx->oldcapo; 846 if (!tjsch->rctx->reverseonestep && tjsch->rctx->stepsleft > 0) { 847 tjsch->rctx->stepsleft--; 848 PetscPrintf(PETSC_COMM_WORLD,"\x1B[35mSkip the step from %D to %D (already checkpointed)\033[0m\n",tjsch->rctx->oldcapo,tjsch->rctx->oldcapo+1); 849 } 850 } 851 if (stack->solution_only || (!stack->solution_only && e->stepnum < stepnum)) { 852 tjsch->recompute = PETSC_TRUE; 853 ierr = ReCompute(ts,tjsch,e->stepnum,stepnum);CHKERRQ(ierr); 854 } 855 if (tjsch->rctx->reverseonestep) tjsch->rctx->reverseonestep = PETSC_FALSE; 856 PetscFunctionReturn(0); 857 } 858 859 #undef __FUNCT__ 860 #define __FUNCT__ "SetTrajTLR" 861 static PetscErrorCode SetTrajTLR(TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X) 862 { 863 Stack *stack = &tjsch->stack; 864 PetscInt store,localstepnum,laststridesize; 865 StackElement e; 866 RevolveCTX *rctx = tjsch->rctx; 867 PetscBool resetrevolve = PETSC_FALSE; 868 PetscErrorCode ierr; 869 870 PetscFunctionBegin; 871 if (!stack->solution_only && stepnum == 0) PetscFunctionReturn(0); 872 if (stack->solution_only && stepnum == tjsch->total_steps) PetscFunctionReturn(0); 873 874 localstepnum = stepnum%tjsch->stride; 875 laststridesize = tjsch->total_steps%tjsch->stride; 876 if (!laststridesize) laststridesize = tjsch->stride; 877 878 if (!tjsch->recompute) { /* skip last stride */ 879 TopLevelStoreAll(ts,tjsch,stepnum,localstepnum,laststridesize,&resetrevolve); 880 if (resetrevolve) PetscPrintf(PETSC_COMM_WORLD,"\x1B[33mDump stack to file\033[0m\n"); 881 else PetscPrintf(PETSC_COMM_WORLD,"\x1B[33mDump a single point to file\033[0m\n"); 882 /* different starting points for last stride between solutin_only and !solutin_only */ 883 if (!stack->solution_only && !tjsch->save_stack && stepnum <= tjsch->total_steps-laststridesize) PetscFunctionReturn(0); 884 if (stack->solution_only && !tjsch->save_stack && stepnum < tjsch->total_steps-laststridesize) PetscFunctionReturn(0); 885 } 886 887 if (tjsch->save_stack && resetrevolve) { 888 revolve_reset(); 889 revolve_create_offline(tjsch->stride,tjsch->max_cps_ram); 890 rctx = tjsch->rctx; 891 rctx->check = 0; 892 rctx->capo = 0; 893 rctx->fine = tjsch->stride; 894 tjsch->rctx->reverseonestep = PETSC_FALSE; 895 PetscFunctionReturn(0); 896 } 897 ierr = ApplyRevolve(tjsch->stype,tjsch->rctx,tjsch->total_steps,stepnum,localstepnum,&store);CHKERRQ(ierr); 898 if (store == 1) { 899 if (localstepnum < stack->top) SETERRQ(tjsch->comm,PETSC_ERR_MEMC,"Illegal modification of a non-top stack element"); 900 ierr = ElementCreate(ts,stack,&e,stepnum,time,X);CHKERRQ(ierr); 901 ierr = StackPush(stack,e);CHKERRQ(ierr); 902 } 903 904 PetscFunctionReturn(0); 905 } 906 907 #undef __FUNCT__ 908 #define __FUNCT__ "GetTrajTLR" 909 static PetscErrorCode GetTrajTLR(TS ts,TJScheduler *tjsch,PetscInt stepnum) 910 { 911 Stack *stack = &tjsch->stack; 912 PetscInt whattodo,shift; 913 PetscInt localstepnum,stridenum,laststridesize,store; 914 PetscReal stepsize; 915 StackElement e; 916 PetscErrorCode ierr; 917 918 PetscFunctionBegin; 919 localstepnum = stepnum%tjsch->stride; 920 stridenum = stepnum/tjsch->stride; 921 if (stepnum == tjsch->total_steps) { 922 ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr); 923 ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr); 924 if ( tjsch->rctx->reverseonestep) tjsch->rctx->reverseonestep = PETSC_FALSE; 925 PetscFunctionReturn(0); 926 } 927 laststridesize = tjsch->total_steps%tjsch->stride; 928 if (!laststridesize) laststridesize = tjsch->stride; 929 if (stack->solution_only) { 930 /* fill stack */ 931 if (localstepnum == 0 && stepnum <= tjsch->total_steps-laststridesize) { 932 if (tjsch->save_stack) { 933 PetscPrintf(PETSC_COMM_WORLD,"\x1B[33mLoad stack from file\033[0m\n"); 934 ierr = StackLoadAll(ts,stack,stridenum);CHKERRQ(ierr); 935 revolve_reset(); 936 revolve_create_offline(tjsch->stride,tjsch->max_cps_ram); 937 tjsch->rctx->check = 0; 938 tjsch->rctx->capo = 0; 939 tjsch->rctx->fine = tjsch->stride; 940 whattodo = 0; 941 while(whattodo!=3) { /* stupid revolve */ 942 whattodo = revolve_action(&tjsch->rctx->check,&tjsch->rctx->capo,&tjsch->rctx->fine,tjsch->rctx->snaps_in,&tjsch->rctx->info,&tjsch->rctx->where); 943 } 944 tjsch->recompute = PETSC_TRUE; 945 tjsch->skip_trajectory = PETSC_TRUE; 946 ierr = ReCompute(ts,tjsch,stridenum*tjsch->stride-1,stridenum*tjsch->stride);CHKERRQ(ierr); 947 tjsch->skip_trajectory = PETSC_FALSE; 948 } else { 949 PetscPrintf(PETSC_COMM_WORLD,"\x1B[33mLoad a single point from file\033[0m\n"); 950 ierr = LoadSingle(ts,stack,stridenum);CHKERRQ(ierr); 951 revolve_reset(); 952 revolve_create_offline(tjsch->stride,tjsch->max_cps_ram); 953 tjsch->rctx->check = 0; 954 tjsch->rctx->capo = 0; 955 tjsch->rctx->fine = tjsch->stride; 956 tjsch->recompute = PETSC_TRUE; 957 ierr = ReCompute(ts,tjsch,(stridenum-1)*tjsch->stride,stridenum*tjsch->stride);CHKERRQ(ierr); 958 } 959 PetscFunctionReturn(0); 960 } 961 /* restore a checkpoint */ 962 ierr = StackTop(stack,&e);CHKERRQ(ierr); 963 ierr = UpdateTS(ts,stack,e);CHKERRQ(ierr); 964 /* start with restoring a checkpoint */ 965 tjsch->rctx->capo = stepnum; 966 tjsch->rctx->oldcapo = tjsch->rctx->capo; 967 shift = stepnum-localstepnum; 968 whattodo = revolve_action(&tjsch->rctx->check,&tjsch->rctx->capo,&tjsch->rctx->fine,tjsch->rctx->snaps_in,&tjsch->rctx->info,&tjsch->rctx->where); 969 printwhattodo(whattodo,tjsch->rctx,shift); 970 tjsch->recompute = PETSC_TRUE; 971 ierr = ReCompute(ts,tjsch,e->stepnum,stepnum);CHKERRQ(ierr); 972 if (e->stepnum+1 == stepnum) { 973 ierr = StackPop(stack,&e);CHKERRQ(ierr); 974 ierr = ElementDestroy(stack,e);CHKERRQ(ierr); 975 } 976 } else { 977 /* fill stack with info */ 978 if (localstepnum == 0 && tjsch->total_steps-stepnum >= laststridesize) { 979 if (tjsch->save_stack) { 980 PetscPrintf(PETSC_COMM_WORLD,"\x1B[33mLoad stack from file\033[0m\n"); 981 ierr = StackLoadAll(ts,stack,stridenum);CHKERRQ(ierr); 982 revolve_reset(); 983 revolve_create_offline(tjsch->stride,tjsch->max_cps_ram); 984 tjsch->rctx->check = 0; 985 tjsch->rctx->capo = 0; 986 tjsch->rctx->fine = tjsch->stride; 987 whattodo = 0; 988 while(whattodo!=3) { /* stupid revolve */ 989 whattodo = revolve_action(&tjsch->rctx->check,&tjsch->rctx->capo,&tjsch->rctx->fine,tjsch->rctx->snaps_in,&tjsch->rctx->info,&tjsch->rctx->where); 990 } 991 } else { 992 PetscPrintf(PETSC_COMM_WORLD,"\x1B[33mLoad a single point from file\033[0m\n"); 993 ierr = LoadSingle(ts,stack,stridenum);CHKERRQ(ierr); 994 revolve_reset(); 995 revolve_create_offline(tjsch->stride,tjsch->max_cps_ram); 996 tjsch->rctx->check = 0; 997 tjsch->rctx->capo = 0; 998 tjsch->rctx->fine = tjsch->stride; 999 shift = stepnum-localstepnum; 1000 whattodo = revolve_action(&tjsch->rctx->check,&tjsch->rctx->capo,&tjsch->rctx->fine,tjsch->rctx->snaps_in,&tjsch->rctx->info,&tjsch->rctx->where); 1001 printwhattodo(whattodo,tjsch->rctx,shift); 1002 ierr = ElementCreate(ts,stack,&e,(stridenum-1)*tjsch->stride+1,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 1003 ierr = StackPush(stack,e);CHKERRQ(ierr); 1004 tjsch->recompute = PETSC_TRUE; 1005 ierr = ReCompute(ts,tjsch,e->stepnum,stridenum*tjsch->stride);CHKERRQ(ierr); 1006 if ( tjsch->rctx->reverseonestep) { /* ready for the reverse step */ 1007 ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr); 1008 ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr); 1009 tjsch->rctx->reverseonestep = PETSC_FALSE; 1010 } 1011 } 1012 PetscFunctionReturn(0); 1013 } 1014 /* restore a checkpoint */ 1015 ierr = StackTop(stack,&e);CHKERRQ(ierr); 1016 ierr = UpdateTS(ts,stack,e);CHKERRQ(ierr); 1017 /* 2 revolve actions: restore a checkpoint and then advance */ 1018 ierr = ApplyRevolve(tjsch->stype,tjsch->rctx,tjsch->total_steps,stepnum,localstepnum,&store);CHKERRQ(ierr); 1019 if (!tjsch->rctx->reverseonestep && tjsch->rctx->stepsleft > 0) { 1020 tjsch->rctx->stepsleft--; 1021 PetscPrintf(PETSC_COMM_WORLD,"\x1B[35mSkip the step from %D to %D (already checkpointed)\033[0m\n",stepnum-localstepnum+tjsch->rctx->oldcapo,stepnum-localstepnum+tjsch->rctx->oldcapo+1); 1022 } 1023 if (e->stepnum < stepnum) { 1024 tjsch->recompute = PETSC_TRUE; 1025 ierr = ReCompute(ts,tjsch,e->stepnum,stepnum);CHKERRQ(ierr); 1026 } 1027 if (e->stepnum == stepnum) { 1028 ierr = StackPop(stack,&e);CHKERRQ(ierr); 1029 ierr = ElementDestroy(stack,e);CHKERRQ(ierr); 1030 } 1031 } 1032 if (tjsch->rctx->reverseonestep) tjsch->rctx->reverseonestep = PETSC_FALSE; 1033 PetscFunctionReturn(0); 1034 } 1035 1036 #undef __FUNCT__ 1037 #define __FUNCT__ "SetTrajTLTR" 1038 static PetscErrorCode SetTrajTLTR(TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X) 1039 { 1040 } 1041 1042 #undef __FUNCT__ 1043 #define __FUNCT__ "GetTrajTLTR" 1044 static PetscErrorCode GetTrajTLTR(TS ts,TJScheduler *tjsch,PetscInt stepnum) 1045 { 1046 } 1047 1048 #undef __FUNCT__ 1049 #define __FUNCT__ "SetTrajRMS" 1050 static PetscErrorCode SetTrajRMS(TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X) 1051 { 1052 Stack *stack = &tjsch->stack; 1053 PetscInt store; 1054 StackElement e; 1055 PetscErrorCode ierr; 1056 1057 PetscFunctionBegin; 1058 if (!stack->solution_only && stepnum == 0) PetscFunctionReturn(0); 1059 if (stack->solution_only && stepnum == tjsch->total_steps) PetscFunctionReturn(0); 1060 ierr = ApplyRevolve(tjsch->stype,tjsch->rctx,tjsch->total_steps,stepnum,stepnum,&store);CHKERRQ(ierr); 1061 if (store == 1){ 1062 if (stepnum < stack->top) SETERRQ(tjsch->comm,PETSC_ERR_MEMC,"Illegal modification of a non-top stack element"); 1063 ierr = ElementCreate(ts,stack,&e,stepnum,time,X);CHKERRQ(ierr); 1064 ierr = StackPush(stack,e);CHKERRQ(ierr); 1065 } else if (store == 2) { 1066 ierr = DumpSingle(ts,stack,tjsch->rctx->check);CHKERRQ(ierr); 1067 } 1068 PetscFunctionReturn(0); 1069 } 1070 1071 #undef __FUNCT__ 1072 #define __FUNCT__ "GetTrajRMS" 1073 static PetscErrorCode GetTrajRMS(TS ts,TJScheduler *tjsch,PetscInt stepnum) 1074 { 1075 Stack *stack = &tjsch->stack; 1076 PetscInt whattodo,shift; 1077 PetscInt restart; 1078 PetscBool ondisk; 1079 PetscReal stepsize; 1080 StackElement e; 1081 PetscErrorCode ierr; 1082 1083 PetscFunctionBegin; 1084 if (stepnum == 0 || stepnum == tjsch->total_steps) { 1085 ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr); 1086 ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr); 1087 if (tjsch->rctx->reverseonestep) tjsch->rctx->reverseonestep = PETSC_FALSE; 1088 PetscFunctionReturn(0); 1089 } 1090 tjsch->rctx->capo = stepnum; 1091 tjsch->rctx->oldcapo = tjsch->rctx->capo; 1092 shift = 0; 1093 whattodo = revolve_action(&tjsch->rctx->check,&tjsch->rctx->capo,&tjsch->rctx->fine,tjsch->rctx->snaps_in,&tjsch->rctx->info,&tjsch->rctx->where); /* whattodo=restore */ 1094 printwhattodo(whattodo,tjsch->rctx,shift); 1095 /* restore a checkpoint */ 1096 restart = tjsch->rctx->capo; 1097 if (!tjsch->rctx->where) { 1098 ondisk = PETSC_TRUE; 1099 ierr = LoadSingle(ts,stack,tjsch->rctx->check);CHKERRQ(ierr); 1100 ierr = TSSetTimeStep(ts,ts->ptime_prev-ts->ptime);CHKERRQ(ierr); 1101 } else { 1102 ondisk = PETSC_FALSE; 1103 ierr = StackTop(stack,&e);CHKERRQ(ierr); 1104 ierr = UpdateTS(ts,stack,e);CHKERRQ(ierr); 1105 } 1106 if (!stack->solution_only) { /* whattodo must be 5 or 8 */ 1107 /* ask Revolve what to do next */ 1108 tjsch->rctx->oldcapo = tjsch->rctx->capo; 1109 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*/ 1110 printwhattodo(whattodo,tjsch->rctx,shift); 1111 if (whattodo == 3 || whattodo == 4) tjsch->rctx->reverseonestep = PETSC_TRUE; 1112 if (whattodo == 1) tjsch->rctx->stepsleft = tjsch->rctx->capo-tjsch->rctx->oldcapo; 1113 if (!tjsch->rctx->reverseonestep && tjsch->rctx->stepsleft > 0) { 1114 tjsch->rctx->stepsleft--; 1115 PetscPrintf(PETSC_COMM_WORLD,"\x1B[35mSkip the step from %D to %D (already checkpointed)\033[0m\n",tjsch->rctx->oldcapo,tjsch->rctx->oldcapo+1); 1116 } 1117 restart++; /* skip one step */ 1118 } 1119 if (stack->solution_only || (!stack->solution_only && restart < stepnum)) { 1120 tjsch->recompute = PETSC_TRUE; 1121 ierr = ReCompute(ts,tjsch,restart,stepnum);CHKERRQ(ierr); 1122 } 1123 if (!ondisk && ( (stack->solution_only && e->stepnum+1 == stepnum) || (!stack->solution_only && e->stepnum == stepnum) )) { 1124 ierr = StackPop(stack,&e);CHKERRQ(ierr); 1125 ierr = ElementDestroy(stack,e);CHKERRQ(ierr); 1126 } 1127 if (tjsch->rctx->reverseonestep) tjsch->rctx->reverseonestep = PETSC_FALSE; 1128 PetscFunctionReturn(0); 1129 } 1130 #endif 1131 1132 #undef __FUNCT__ 1133 #define __FUNCT__ "TSTrajectorySet_Memory" 1134 PetscErrorCode TSTrajectorySet_Memory(TSTrajectory tj,TS ts,PetscInt stepnum,PetscReal time,Vec X) 1135 { 1136 TJScheduler *tjsch = (TJScheduler*)tj->data; 1137 PetscErrorCode ierr; 1138 1139 PetscFunctionBegin; 1140 if (!tjsch->recompute) { /* use global stepnum in the forward sweep */ 1141 ierr = TSGetTotalSteps(ts,&stepnum);CHKERRQ(ierr); 1142 } 1143 /* for consistency */ 1144 if (!tjsch->recompute && stepnum == 0) ts->ptime_prev = ts->ptime-ts->time_step; 1145 switch (tjsch->stype) { 1146 case NONE: 1147 ierr = SetTrajN(ts,tjsch,stepnum,time,X);CHKERRQ(ierr); 1148 break; 1149 case TWO_LEVEL_NOREVOLVE: 1150 ierr = SetTrajTLNR(ts,tjsch,stepnum,time,X);CHKERRQ(ierr); 1151 break; 1152 #ifdef PETSC_HAVE_REVOLVE 1153 case TWO_LEVEL_REVOLVE: 1154 ierr = SetTrajTLR(ts,tjsch,stepnum,time,X);CHKERRQ(ierr); 1155 break; 1156 case TWO_LEVEL_TWO_REVOLVE: 1157 ierr = SetTrajTLTR(ts,tjsch,stepnum,time,X);CHKERRQ(ierr); 1158 break; 1159 case REVOLVE_OFFLINE: 1160 ierr = SetTrajROF(ts,tjsch,stepnum,time,X);CHKERRQ(ierr); 1161 break; 1162 case REVOLVE_ONLINE: 1163 ierr = SetTrajRON(ts,tjsch,stepnum,time,X);CHKERRQ(ierr); 1164 break; 1165 case REVOLVE_MULTISTAGE: 1166 ierr = SetTrajRMS(ts,tjsch,stepnum,time,X);CHKERRQ(ierr); 1167 break; 1168 #endif 1169 default: 1170 break; 1171 } 1172 PetscFunctionReturn(0); 1173 } 1174 1175 #undef __FUNCT__ 1176 #define __FUNCT__ "TSTrajectoryGet_Memory" 1177 PetscErrorCode TSTrajectoryGet_Memory(TSTrajectory tj,TS ts,PetscInt stepnum,PetscReal *t) 1178 { 1179 TJScheduler *tjsch = (TJScheduler*)tj->data; 1180 PetscErrorCode ierr; 1181 1182 PetscFunctionBegin; 1183 ierr = TSGetTotalSteps(ts,&stepnum);CHKERRQ(ierr); 1184 if (stepnum == 0) PetscFunctionReturn(0); 1185 switch (tjsch->stype) { 1186 case NONE: 1187 ierr = GetTrajN(ts,tjsch,stepnum);CHKERRQ(ierr); 1188 break; 1189 case TWO_LEVEL_NOREVOLVE: 1190 ierr = GetTrajTLNR(ts,tjsch,stepnum);CHKERRQ(ierr); 1191 break; 1192 #ifdef PETSC_HAVE_REVOLVE 1193 case TWO_LEVEL_REVOLVE: 1194 ierr = GetTrajTLR(ts,tjsch,stepnum);CHKERRQ(ierr); 1195 break; 1196 case TWO_LEVEL_TWO_REVOLVE: 1197 ierr = GetTrajTLTR(ts,tjsch,stepnum);CHKERRQ(ierr); 1198 break; 1199 case REVOLVE_OFFLINE: 1200 ierr = GetTrajROF(ts,tjsch,stepnum);CHKERRQ(ierr); 1201 break; 1202 case REVOLVE_ONLINE: 1203 ierr = GetTrajRON(ts,tjsch,stepnum);CHKERRQ(ierr); 1204 break; 1205 case REVOLVE_MULTISTAGE: 1206 ierr = GetTrajRMS(ts,tjsch,stepnum);CHKERRQ(ierr); 1207 break; 1208 #endif 1209 default: 1210 break; 1211 } 1212 PetscFunctionReturn(0); 1213 } 1214 1215 #undef __FUNCT__ 1216 #define __FUNCT__ "TSTrajectorySetStride_Memory" 1217 PetscErrorCode TSTrajectorySetStride_Memory(TSTrajectory tj,TS ts,PetscInt stride) 1218 { 1219 TJScheduler *tjsch = (TJScheduler*)tj->data; 1220 1221 PetscFunctionBegin; 1222 tjsch->stride = stride; 1223 PetscFunctionReturn(0); 1224 } 1225 1226 #undef __FUNCT__ 1227 #define __FUNCT__ "TSTrajectorySetMaxCpsRAM_Memory" 1228 PetscErrorCode TSTrajectorySetMaxCpsRAM_Memory(TSTrajectory tj,TS ts,PetscInt max_cps_ram) 1229 { 1230 TJScheduler *tjsch = (TJScheduler*)tj->data; 1231 1232 PetscFunctionBegin; 1233 tjsch->max_cps_ram = max_cps_ram; 1234 PetscFunctionReturn(0); 1235 } 1236 1237 #undef __FUNCT__ 1238 #define __FUNCT__ "TSTrajectorySetMaxCpsDisk_Memory" 1239 PetscErrorCode TSTrajectorySetMaxCpsDisk_Memory(TSTrajectory tj,TS ts,PetscInt max_cps_disk) 1240 { 1241 TJScheduler *tjsch = (TJScheduler*)tj->data; 1242 1243 PetscFunctionBegin; 1244 tjsch->max_cps_disk = max_cps_disk; 1245 PetscFunctionReturn(0); 1246 } 1247 1248 #ifdef PETSC_HAVE_REVOLVE 1249 #undef __FUNCT__ 1250 #define __FUNCT__ "TSTrajectorySetRevolveOnline" 1251 PetscErrorCode TSTrajectorySetRevolveOnline(TSTrajectory tj,PetscBool use_online) 1252 { 1253 TJScheduler *tjsch = (TJScheduler*)tj->data; 1254 1255 PetscFunctionBegin; 1256 tjsch->use_online = use_online; 1257 PetscFunctionReturn(0); 1258 } 1259 #endif 1260 1261 #undef __FUNCT__ 1262 #define __FUNCT__ "TSTrajectorySetSaveStack" 1263 PetscErrorCode TSTrajectorySetSaveStack(TSTrajectory tj,PetscBool save_stack) 1264 { 1265 TJScheduler *tjsch = (TJScheduler*)tj->data; 1266 1267 PetscFunctionBegin; 1268 tjsch->save_stack = save_stack; 1269 PetscFunctionReturn(0); 1270 } 1271 1272 #undef __FUNCT__ 1273 #define __FUNCT__ "TSTrajectorySetSolutionOnly" 1274 PetscErrorCode TSTrajectorySetSolutionOnly(TSTrajectory tj,PetscBool solution_only) 1275 { 1276 TJScheduler *tjsch = (TJScheduler*)tj->data; 1277 Stack *stack = &tjsch->stack; 1278 1279 PetscFunctionBegin; 1280 stack->solution_only = solution_only; 1281 PetscFunctionReturn(0); 1282 } 1283 1284 #undef __FUNCT__ 1285 #define __FUNCT__ "TSTrajectorySetFromOptions_Memory" 1286 PetscErrorCode TSTrajectorySetFromOptions_Memory(PetscOptions *PetscOptionsObject,TSTrajectory tj) 1287 { 1288 TJScheduler *tjsch = (TJScheduler*)tj->data; 1289 PetscErrorCode ierr; 1290 1291 PetscFunctionBegin; 1292 ierr = PetscOptionsHead(PetscOptionsObject,"Memory based TS trajectory options");CHKERRQ(ierr); 1293 { 1294 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); 1295 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); 1296 ierr = PetscOptionsInt("-tstrajectory_stride","Stride to save checkpoints to file","TSTrajectorySetStride_Memory",tjsch->stride,&tjsch->stride,NULL);CHKERRQ(ierr); 1297 #ifdef PETSC_HAVE_REVOLVE 1298 ierr = PetscOptionsBool("-tstrajectory_revolve_online","Trick TS trajectory into using online mode of revolve","TSTrajectorySetRevolveOnline",tjsch->use_online,&tjsch->use_online,NULL);CHKERRQ(ierr); 1299 #endif 1300 ierr = PetscOptionsBool("-tstrajectory_save_stack","Save all stack to disk","TSTrajectorySetSaveStack",tjsch->save_stack,&tjsch->save_stack,NULL);CHKERRQ(ierr); 1301 ierr = PetscOptionsBool("-tstrajectory_solution_only","Checkpoint solution only","TSTrajectorySetSolutionOnly",tjsch->stack.solution_only,&tjsch->stack.solution_only,NULL);CHKERRQ(ierr); 1302 } 1303 ierr = PetscOptionsTail();CHKERRQ(ierr); 1304 PetscFunctionReturn(0); 1305 } 1306 1307 #undef __FUNCT__ 1308 #define __FUNCT__ "TSTrajectorySetUp_Memory" 1309 PetscErrorCode TSTrajectorySetUp_Memory(TSTrajectory tj,TS ts) 1310 { 1311 TJScheduler *tjsch = (TJScheduler*)tj->data; 1312 Stack *stack = &tjsch->stack; 1313 #ifdef PETSC_HAVE_REVOLVE 1314 RevolveCTX *rctx; 1315 #endif 1316 PetscInt numY,ns; 1317 PetscBool flg; 1318 PetscErrorCode ierr; 1319 1320 PetscFunctionBegin; 1321 PetscStrcmp(((PetscObject)ts->adapt)->type_name,TSADAPTNONE,&flg); 1322 if (flg) { /* fixed time step */ 1323 tjsch->total_steps = PetscMin(ts->max_steps,(PetscInt)(ceil(ts->max_time/ts->time_step))); 1324 } 1325 if (tjsch->max_cps_ram > 1) stack->stacksize = tjsch->max_cps_ram; 1326 if (tjsch->stride > 1) { /* two level mode works for both fixed time step and adaptive time step */ 1327 if (tjsch->max_cps_disk <= 1 && tjsch->max_cps_ram > 1 && tjsch->max_cps_ram < tjsch->stride-1) { /* use revolve_offline for each stride */ 1328 tjsch->stype = TWO_LEVEL_REVOLVE; 1329 } 1330 if (tjsch->max_cps_disk > 1 && tjsch->max_cps_ram > 1 && tjsch->max_cps_ram <= tjsch->stride-1) { /* use revolve_offline for each stride */ 1331 tjsch->stype = TWO_LEVEL_TWO_REVOLVE; 1332 tjsch->diskstack.stacksize = tjsch->max_cps_disk; 1333 } 1334 if (tjsch->max_cps_disk <= 1 && (tjsch->max_cps_ram >= tjsch->stride || tjsch->max_cps_ram == -1)) { /* checkpoint all for each stride */ 1335 tjsch->stype = TWO_LEVEL_NOREVOLVE; /* can also be handled by TWO_LEVEL_REVOLVE */ 1336 stack->stacksize = tjsch->stride-1; /* need tjsch->stride-1 at most */ 1337 } 1338 } else { 1339 if (flg) { /* fixed time step */ 1340 if (tjsch->max_cps_ram >= tjsch->total_steps-1 || tjsch->max_cps_ram < 1) { /* checkpoint all */ 1341 tjsch->stype = NONE; 1342 stack->stacksize = stack->solution_only ? tjsch->total_steps : tjsch->total_steps-1; 1343 } else { 1344 if (tjsch->max_cps_disk > 1) { /* disk can be used */ 1345 tjsch->stype = REVOLVE_MULTISTAGE; 1346 } else { /* memory only */ 1347 tjsch->stype = REVOLVE_OFFLINE; 1348 } 1349 } 1350 } else { /* adaptive time step */ 1351 tjsch->stype = REVOLVE_ONLINE; 1352 } 1353 #ifdef PETSC_HAVE_REVOLVE 1354 if (tjsch->use_online) { /* trick into online */ 1355 tjsch->stype = REVOLVE_ONLINE; 1356 stack->stacksize = tjsch->max_cps_ram; 1357 } 1358 #endif 1359 } 1360 1361 if (tjsch->stype > TWO_LEVEL_NOREVOLVE) { 1362 #ifndef PETSC_HAVE_REVOLVE 1363 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."); 1364 #else 1365 if (tjsch->stype == TWO_LEVEL_REVOLVE) revolve_create_offline(tjsch->stride,tjsch->max_cps_ram); 1366 else if (tjsch->stype == TWO_LEVEL_TWO_REVOLVE) { 1367 revolve_create_offline(tjsch->stride,tjsch->max_cps_ram); 1368 ns = (tjsch->total_steps+tjsch->stride-1)/tjsch->stride; 1369 revolve2_create_offline(ns,tjsch->max_cps_disk); 1370 } else if (tjsch->stype == REVOLVE_OFFLINE) revolve_create_offline(tjsch->total_steps,tjsch->max_cps_ram); 1371 else if (tjsch->stype == REVOLVE_ONLINE) revolve_create_online(tjsch->max_cps_ram); 1372 else if (tjsch->stype == REVOLVE_MULTISTAGE) revolve_create_multistage(tjsch->total_steps,tjsch->max_cps_ram+tjsch->max_cps_disk,tjsch->max_cps_ram); 1373 1374 ierr = PetscCalloc1(1,&rctx);CHKERRQ(ierr); 1375 rctx->snaps_in = tjsch->max_cps_ram; /* for theta methods snaps_in=2*max_cps_ram */ 1376 rctx->reverseonestep = PETSC_FALSE; 1377 rctx->check = 0; 1378 rctx->oldcapo = 0; 1379 rctx->capo = 0; 1380 rctx->info = 2; 1381 rctx->fine = (tjsch->stride > 1) ? tjsch->stride : tjsch->total_steps; 1382 1383 tjsch->rctx = rctx; 1384 if (tjsch->stype == REVOLVE_ONLINE) rctx->fine = -1; 1385 #endif 1386 } 1387 1388 tjsch->recompute = PETSC_FALSE; 1389 tjsch->comm = PetscObjectComm((PetscObject)ts); 1390 ierr = TSGetStages(ts,&numY,PETSC_IGNORE);CHKERRQ(ierr); 1391 ierr = StackCreate(stack,stack->stacksize,numY);CHKERRQ(ierr); 1392 PetscFunctionReturn(0); 1393 } 1394 1395 #undef __FUNCT__ 1396 #define __FUNCT__ "TSTrajectoryDestroy_Memory" 1397 PETSC_EXTERN PetscErrorCode TSTrajectoryDestroy_Memory(TSTrajectory tj) 1398 { 1399 TJScheduler *tjsch = (TJScheduler*)tj->data; 1400 PetscErrorCode ierr; 1401 1402 PetscFunctionBegin; 1403 if (tjsch->stype > TWO_LEVEL_NOREVOLVE) { 1404 #ifdef PETSC_HAVE_REVOLVE 1405 revolve_reset(); 1406 if (tjsch->stype == TWO_LEVEL_TWO_REVOLVE) revolve2_reset(); 1407 #endif 1408 } 1409 ierr = StackDestroy(&tjsch->stack);CHKERRQ(ierr); 1410 #ifdef PETSC_HAVE_REVOLVE 1411 if (tjsch->stype > TWO_LEVEL_NOREVOLVE) { 1412 ierr = PetscFree(tjsch->rctx);CHKERRQ(ierr); 1413 } 1414 #endif 1415 ierr = PetscFree(tjsch);CHKERRQ(ierr); 1416 PetscFunctionReturn(0); 1417 } 1418 1419 /*MC 1420 TSTRAJECTORYMEMORY - Stores each solution of the ODE/ADE in memory 1421 1422 Level: intermediate 1423 1424 .seealso: TSTrajectoryCreate(), TS, TSTrajectorySetType() 1425 1426 M*/ 1427 #undef __FUNCT__ 1428 #define __FUNCT__ "TSTrajectoryCreate_Memory" 1429 PETSC_EXTERN PetscErrorCode TSTrajectoryCreate_Memory(TSTrajectory tj,TS ts) 1430 { 1431 TJScheduler *tjsch; 1432 PetscErrorCode ierr; 1433 1434 PetscFunctionBegin; 1435 tj->ops->set = TSTrajectorySet_Memory; 1436 tj->ops->get = TSTrajectoryGet_Memory; 1437 tj->ops->setup = TSTrajectorySetUp_Memory; 1438 tj->ops->destroy = TSTrajectoryDestroy_Memory; 1439 tj->ops->setfromoptions = TSTrajectorySetFromOptions_Memory; 1440 1441 ierr = PetscCalloc1(1,&tjsch);CHKERRQ(ierr); 1442 tjsch->stype = NONE; 1443 tjsch->max_cps_ram = -1; /* -1 indicates that it is not set */ 1444 tjsch->max_cps_disk = -1; /* -1 indicates that it is not set */ 1445 tjsch->stride = 0; /* if not zero, two-level checkpointing will be used */ 1446 #ifdef PETSC_HAVE_REVOLVE 1447 tjsch->use_online = PETSC_FALSE; 1448 #endif 1449 tjsch->save_stack = PETSC_TRUE; 1450 1451 tjsch->stack.solution_only = PETSC_TRUE; 1452 1453 tj->data = tjsch; 1454 1455 PetscFunctionReturn(0); 1456 } 1457