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__ "TopLevelStore" 429 static PetscErrorCode TopLevelStore(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 /* make sure saved checkpoint id starts from 1 438 skip last stride when using stridenum+1 439 skip first stride when using stridenum */ 440 if (stack->solution_only) { 441 if (tjsch->save_stack) { 442 if (localstepnum == tjsch->stride-1 && stepnum < tjsch->total_steps-laststridesize) { /* current step will be saved without going through stack */ 443 ierr = StackDumpAll(ts,stack,stridenum+1);CHKERRQ(ierr); 444 PetscPrintf(PETSC_COMM_WORLD,"\x1B[33mDump stack to file\033[0m\n"); 445 if (reset) *reset = PETSC_TRUE; 446 } 447 } else { 448 if (!tjsch->recompute && localstepnum == 0 && stepnum < tjsch->total_steps-laststridesize) { 449 ierr = DumpSingle(ts,stack,stridenum+1);CHKERRQ(ierr); 450 PetscPrintf(PETSC_COMM_WORLD,"\x1B[33mDump a single point (solution) to file\033[0m\n"); 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 PetscPrintf(PETSC_COMM_WORLD,"\x1B[33mDump stack to file\033[0m\n"); 458 if (reset) *reset = PETSC_TRUE; 459 } 460 } else { 461 if (!tjsch->recompute && localstepnum == 1 && stepnum < tjsch->total_steps-laststridesize) { 462 ierr = DumpSingle(ts,stack,stridenum+1);CHKERRQ(ierr); 463 PetscPrintf(PETSC_COMM_WORLD,"\x1B[33mDump a single point (solution+stages) to file\033[0m\n"); 464 } 465 } 466 } 467 PetscFunctionReturn(0); 468 } 469 470 #undef __FUNCT__ 471 #define __FUNCT__ "SetTrajN" 472 static PetscErrorCode SetTrajN(TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X) 473 { 474 Stack *stack = &tjsch->stack; 475 StackElement e; 476 PetscErrorCode ierr; 477 478 PetscFunctionBegin; 479 /* skip the last two steps of each stride or the whole interval */ 480 if (stack->solution_only && (stepnum >= tjsch->total_steps-1 || tjsch->recompute)) PetscFunctionReturn(0); //? 481 /* skip the first and the last steps of each stride or the whole interval */ 482 if (!stack->solution_only && (stepnum == 0 || stepnum == tjsch->total_steps)) PetscFunctionReturn(0); 483 484 if (stepnum < stack->top) SETERRQ(tjsch->comm,PETSC_ERR_MEMC,"Illegal modification of a non-top stack element"); 485 ierr = ElementCreate(ts,stack,&e,stepnum,time,X);CHKERRQ(ierr); 486 ierr = StackPush(stack,e);CHKERRQ(ierr); 487 PetscFunctionReturn(0); 488 } 489 490 #undef __FUNCT__ 491 #define __FUNCT__ "GetTrajN" 492 static PetscErrorCode GetTrajN(TS ts,TJScheduler *tjsch,PetscInt stepnum) 493 { 494 Stack *stack = &tjsch->stack; 495 PetscReal stepsize; 496 StackElement e; 497 PetscErrorCode ierr; 498 499 PetscFunctionBegin; 500 if (stepnum == tjsch->total_steps) { 501 ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr); 502 ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr); 503 PetscFunctionReturn(0); 504 } 505 /* restore a checkpoint */ 506 ierr = StackTop(stack,&e);CHKERRQ(ierr); 507 ierr = UpdateTS(ts,stack,e);CHKERRQ(ierr); 508 if (stack->solution_only) {/* recompute one step */ 509 tjsch->recompute = PETSC_TRUE; 510 ierr = ReCompute(ts,tjsch,e->stepnum,stepnum);CHKERRQ(ierr); 511 } 512 ierr = StackPop(stack,&e);CHKERRQ(ierr); 513 ierr = ElementDestroy(stack,e);CHKERRQ(ierr); 514 PetscFunctionReturn(0); 515 } 516 517 #undef __FUNCT__ 518 #define __FUNCT__ "SetTrajTLNR" 519 static PetscErrorCode SetTrajTLNR(TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X) 520 { 521 Stack *stack = &tjsch->stack; 522 PetscInt localstepnum,laststridesize; 523 StackElement e; 524 PetscErrorCode ierr; 525 526 PetscFunctionBegin; 527 if (!stack->solution_only && stepnum == 0) PetscFunctionReturn(0); 528 if (stack->solution_only && stepnum == tjsch->total_steps) PetscFunctionReturn(0); 529 if (tjsch->save_stack && tjsch->recompute) PetscFunctionReturn(0); 530 531 localstepnum = stepnum%tjsch->stride; 532 /* (stride size-1) checkpoints are saved in each stride; an extra point is added by StackDumpAll() */ 533 laststridesize = tjsch->total_steps%tjsch->stride; 534 if (!laststridesize) laststridesize = tjsch->stride; 535 536 if (!tjsch->recompute) { 537 ierr = TopLevelStore(ts,tjsch,stepnum,localstepnum,laststridesize,PETSC_NULL);CHKERRQ(ierr); 538 if (!tjsch->save_stack && stepnum < tjsch->total_steps-laststridesize) PetscFunctionReturn(0); 539 } 540 if (!stack->solution_only && localstepnum == 0) PetscFunctionReturn(0); /* skip last point in each stride at recompute stage or last stride */ 541 if (stack->solution_only && localstepnum == tjsch->stride-1) PetscFunctionReturn(0); /* skip last step in each stride at recompute stage or last stride */ 542 543 ierr = ElementCreate(ts,stack,&e,stepnum,time,X);CHKERRQ(ierr); 544 ierr = StackPush(stack,e);CHKERRQ(ierr); 545 PetscFunctionReturn(0); 546 } 547 548 #undef __FUNCT__ 549 #define __FUNCT__ "GetTrajTLNR" 550 static PetscErrorCode GetTrajTLNR(TS ts,TJScheduler *tjsch,PetscInt stepnum) 551 { 552 Stack *stack = &tjsch->stack; 553 PetscInt id,localstepnum,laststridesize; 554 PetscReal stepsize; 555 StackElement e; 556 PetscErrorCode ierr; 557 558 PetscFunctionBegin; 559 if (stepnum == tjsch->total_steps) { 560 ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr); 561 ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr); 562 PetscFunctionReturn(0); 563 } 564 565 localstepnum = stepnum%tjsch->stride; 566 laststridesize = tjsch->total_steps%tjsch->stride; 567 if (!laststridesize) laststridesize = tjsch->stride; 568 if (stack->solution_only) { 569 /* fill stack with info */ 570 if (localstepnum == 0 && tjsch->total_steps-stepnum >= laststridesize) { 571 id = stepnum/tjsch->stride; 572 if (tjsch->save_stack) { 573 ierr = StackLoadAll(ts,stack,id);CHKERRQ(ierr); 574 tjsch->recompute = PETSC_TRUE; 575 tjsch->skip_trajectory = PETSC_TRUE; 576 ierr = ReCompute(ts,tjsch,id*tjsch->stride-1,id*tjsch->stride);CHKERRQ(ierr); 577 tjsch->skip_trajectory = PETSC_FALSE; 578 } else { 579 ierr = LoadSingle(ts,stack,id);CHKERRQ(ierr); 580 tjsch->recompute = PETSC_TRUE; 581 ierr = ReCompute(ts,tjsch,(id-1)*tjsch->stride,id*tjsch->stride);CHKERRQ(ierr); 582 } 583 PetscFunctionReturn(0); 584 } 585 /* restore a checkpoint */ 586 ierr = StackPop(stack,&e);CHKERRQ(ierr); 587 ierr = UpdateTS(ts,stack,e);CHKERRQ(ierr); 588 tjsch->recompute = PETSC_TRUE; 589 tjsch->skip_trajectory = PETSC_TRUE; 590 ierr = ReCompute(ts,tjsch,e->stepnum,stepnum);CHKERRQ(ierr); 591 tjsch->skip_trajectory = PETSC_FALSE; 592 ierr = ElementDestroy(stack,e);CHKERRQ(ierr); 593 } else { 594 /* fill stack with info */ 595 if (localstepnum == 0 && tjsch->total_steps-stepnum >= laststridesize) { 596 id = stepnum/tjsch->stride; 597 if (tjsch->save_stack) { 598 ierr = StackLoadAll(ts,stack,id);CHKERRQ(ierr); 599 } else { 600 ierr = LoadSingle(ts,stack,id);CHKERRQ(ierr); 601 ierr = ElementCreate(ts,stack,&e,(id-1)*tjsch->stride+1,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 602 ierr = StackPush(stack,e);CHKERRQ(ierr); 603 tjsch->recompute = PETSC_TRUE; 604 ierr = ReCompute(ts,tjsch,e->stepnum,id*tjsch->stride);CHKERRQ(ierr); 605 } 606 PetscFunctionReturn(0); 607 } 608 /* restore a checkpoint */ 609 ierr = StackPop(stack,&e);CHKERRQ(ierr); 610 ierr = UpdateTS(ts,stack,e);CHKERRQ(ierr); 611 ierr = ElementDestroy(stack,e);CHKERRQ(ierr); 612 } 613 PetscFunctionReturn(0); 614 } 615 616 #ifdef PETSC_HAVE_REVOLVE 617 static void printwhattodo(PetscInt whattodo,RevolveCTX *rctx,PetscInt shift) 618 { 619 switch(whattodo) { 620 case 1: 621 PetscPrintf(PETSC_COMM_WORLD,"\x1B[35mAdvance from %D to %D\033[0m\n",rctx->oldcapo+shift,rctx->capo+shift); 622 break; 623 case 2: 624 PetscPrintf(PETSC_COMM_WORLD,"\x1B[35mStore in checkpoint number %D (located in RAM)\033[0m\n",rctx->check); 625 break; 626 case 3: 627 PetscPrintf(PETSC_COMM_WORLD,"\x1B[35mFirst turn: Initialize adjoints and reverse first step\033[0m\n"); 628 break; 629 case 4: 630 PetscPrintf(PETSC_COMM_WORLD,"\x1B[35mForward and reverse one step\033[0m\n"); 631 break; 632 case 5: 633 PetscPrintf(PETSC_COMM_WORLD,"\x1B[35mRestore in checkpoint number %D (located in RAM)\033[0m\n",rctx->check); 634 break; 635 case 7: 636 PetscPrintf(PETSC_COMM_WORLD,"\x1B[35mStore in checkpoint number %D (located on disk)\033[0m\n",rctx->check); 637 break; 638 case 8: 639 PetscPrintf(PETSC_COMM_WORLD,"\x1B[35mRestore in checkpoint number %D (located on disk)\033[0m\n",rctx->check); 640 break; 641 case -1: 642 PetscPrintf(PETSC_COMM_WORLD,"\x1B[35mError!"); 643 break; 644 } 645 } 646 647 #undef __FUNCT__ 648 #define __FUNCT__ "ApplyRevolve" 649 static PetscErrorCode ApplyRevolve(SchedulerType stype,RevolveCTX *rctx,PetscInt total_steps,PetscInt stepnum,PetscInt localstepnum,PetscInt *store) 650 { 651 PetscInt shift,whattodo; 652 653 PetscFunctionBegin; 654 *store = 0; 655 // if (rctx->reverseonestep && stepnum == total_steps) { /* intermediate information is ready inside TS, this happens at last time step */ 656 // rctx->reverseonestep = PETSC_FALSE; 657 // PetscFunctionReturn(0); 658 // } 659 if (rctx->stepsleft > 0) { /* advance the solution without checkpointing anything as Revolve requires */ 660 rctx->stepsleft--; 661 PetscFunctionReturn(0); 662 } 663 /* let Revolve determine what to do next */ 664 shift = stepnum-localstepnum; 665 rctx->oldcapo = rctx->capo; 666 rctx->capo = localstepnum; 667 whattodo = revolve_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where); 668 if (stype == REVOLVE_ONLINE && whattodo == 8) whattodo = 5; 669 if (stype == REVOLVE_ONLINE && whattodo == 7) whattodo = 2; 670 printwhattodo(whattodo,rctx,shift); 671 if (whattodo == -1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_LIB,"Error in the Revolve library"); 672 if (whattodo == 1) { /* advance some time steps */ 673 if (stype == REVOLVE_ONLINE && rctx->capo >= total_steps-1) { 674 revolve_turn(total_steps,&rctx->capo,&rctx->fine); 675 printwhattodo(whattodo,rctx,shift); 676 } 677 rctx->stepsleft = rctx->capo-rctx->oldcapo-1; 678 } 679 if (whattodo == 3 || whattodo == 4) { /* ready for a reverse step */ 680 rctx->reverseonestep = PETSC_TRUE; 681 } 682 if (whattodo == 5) { /* restore a checkpoint and ask Revolve what to do next */ 683 rctx->oldcapo = rctx->capo; 684 whattodo = revolve_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where); /* must return 1 or 3 or 4*/ 685 printwhattodo(whattodo,rctx,shift); 686 if (whattodo == 3 || whattodo == 4) rctx->reverseonestep = PETSC_TRUE; 687 if (whattodo == 1) rctx->stepsleft = rctx->capo-rctx->oldcapo; 688 } 689 if (whattodo == 7) { /* save the checkpoint to disk */ 690 *store = 2; 691 rctx->oldcapo = rctx->capo; 692 whattodo = revolve_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where); /* must return 1 */ 693 printwhattodo(whattodo,rctx,shift); 694 rctx->stepsleft = rctx->capo-rctx->oldcapo-1; 695 } 696 if (whattodo == 2) { /* store a checkpoint to RAM and ask Revolve how many time steps to advance next */ 697 *store = 1; 698 rctx->oldcapo = rctx->capo; 699 whattodo = revolve_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where); /* must return 1 */ 700 printwhattodo(whattodo,rctx,shift); 701 if (stype == REVOLVE_ONLINE && rctx->capo >= total_steps-1) { 702 revolve_turn(total_steps,&rctx->capo,&rctx->fine); 703 printwhattodo(whattodo,rctx,shift); 704 } 705 rctx->stepsleft = rctx->capo-rctx->oldcapo-1; 706 } 707 PetscFunctionReturn(0); 708 } 709 710 #undef __FUNCT__ 711 #define __FUNCT__ "TopLevelRevolve" 712 static PetscErrorCode TopLevelRevolve(TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscInt localstepnum,PetscInt laststridesize,PetscBool *reset) 713 { 714 Stack *stack = &tjsch->stack; 715 PetscInt stridenum; 716 PetscErrorCode ierr; 717 718 PetscFunctionBegin; 719 stridenum = stepnum/tjsch->stride; 720 /* make sure saved checkpoint id starts from 1 721 skip last stride when using stridenum+1 722 skip first stride when using stridenum */ 723 if (stack->solution_only) { 724 if (tjsch->save_stack) { 725 if (localstepnum == tjsch->stride-1 && stepnum < tjsch->total_steps-laststridesize) { /* current step will be saved without going through stack */ 726 ierr = StackDumpAll(ts,stack,stridenum+1);CHKERRQ(ierr); 727 PetscPrintf(PETSC_COMM_WORLD,"\x1B[33mDump stack to file\033[0m\n"); 728 if (reset) *reset = PETSC_TRUE; 729 } 730 } else { 731 if (!tjsch->recompute && localstepnum == 0 && stepnum < tjsch->total_steps-laststridesize) { 732 ierr = DumpSingle(ts,stack,stridenum+1);CHKERRQ(ierr); 733 PetscPrintf(PETSC_COMM_WORLD,"\x1B[33mDump a single point (solution) to file\033[0m\n"); 734 } 735 } 736 } else { 737 if (tjsch->save_stack) { 738 if (localstepnum == 0 && stepnum < tjsch->total_steps && stepnum != 0) { 739 ierr = StackDumpAll(ts,stack,stridenum);CHKERRQ(ierr); 740 PetscPrintf(PETSC_COMM_WORLD,"\x1B[33mDump stack to file\033[0m\n"); 741 if (reset) *reset = PETSC_TRUE; 742 } 743 } else { 744 if (!tjsch->recompute && localstepnum == 1 && stepnum < tjsch->total_steps-laststridesize) { 745 ierr = DumpSingle(ts,stack,stridenum+1);CHKERRQ(ierr); 746 PetscPrintf(PETSC_COMM_WORLD,"\x1B[33mDump a single point (solution+stages) to file\033[0m\n"); 747 } 748 } 749 } 750 PetscFunctionReturn(0); 751 } 752 753 #undef __FUNCT__ 754 #define __FUNCT__ "SetTrajROF" 755 static PetscErrorCode SetTrajROF(TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X) 756 { 757 Stack *stack = &tjsch->stack; 758 PetscInt store; 759 StackElement e; 760 PetscErrorCode ierr; 761 762 PetscFunctionBegin; 763 if (!stack->solution_only && stepnum == 0) PetscFunctionReturn(0); 764 if (stack->solution_only && stepnum == tjsch->total_steps) PetscFunctionReturn(0); 765 ierr = ApplyRevolve(tjsch->stype,tjsch->rctx,tjsch->total_steps,stepnum,stepnum,&store);CHKERRQ(ierr); 766 if (store == 1) { 767 if (stepnum < stack->top) SETERRQ(tjsch->comm,PETSC_ERR_MEMC,"Illegal modification of a non-top stack element"); 768 ierr = ElementCreate(ts,stack,&e,stepnum,time,X);CHKERRQ(ierr); 769 ierr = StackPush(stack,e);CHKERRQ(ierr); 770 } 771 PetscFunctionReturn(0); 772 } 773 774 #undef __FUNCT__ 775 #define __FUNCT__ "GetTrajROF" 776 static PetscErrorCode GetTrajROF(TS ts,TJScheduler *tjsch,PetscInt stepnum) 777 { 778 Stack *stack = &tjsch->stack; 779 PetscInt whattodo,shift,store; 780 PetscReal stepsize; 781 StackElement e; 782 PetscErrorCode ierr; 783 784 PetscFunctionBegin; 785 if (stepnum == 0 || stepnum == tjsch->total_steps) { 786 ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr); 787 ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr); 788 if (tjsch->rctx->reverseonestep) tjsch->rctx->reverseonestep = PETSC_FALSE; 789 PetscFunctionReturn(0); 790 } 791 /* restore a checkpoint */ 792 ierr = StackTop(stack,&e);CHKERRQ(ierr); 793 ierr = UpdateTS(ts,stack,e);CHKERRQ(ierr); 794 if (stack->solution_only) { /* start with restoring a checkpoint */ 795 tjsch->rctx->capo = stepnum; 796 tjsch->rctx->oldcapo = tjsch->rctx->capo; 797 shift = 0; 798 whattodo = revolve_action(&tjsch->rctx->check,&tjsch->rctx->capo,&tjsch->rctx->fine,tjsch->rctx->snaps_in,&tjsch->rctx->info,&tjsch->rctx->where); 799 printwhattodo(whattodo,tjsch->rctx,shift); 800 } else { /* 2 revolve actions: restore a checkpoint and then advance */ 801 ierr = ApplyRevolve(tjsch->stype,tjsch->rctx,tjsch->total_steps,stepnum,stepnum,&store);CHKERRQ(ierr); 802 if (!tjsch->rctx->reverseonestep && tjsch->rctx->stepsleft > 0) { 803 tjsch->rctx->stepsleft--; 804 PetscPrintf(PETSC_COMM_WORLD,"\x1B[35mSkip the step from %D to %D (already checkpointed)\033[0m\n",tjsch->rctx->oldcapo,tjsch->rctx->oldcapo+1); 805 } 806 } 807 if (stack->solution_only || (!stack->solution_only && e->stepnum < stepnum)) { 808 tjsch->recompute = PETSC_TRUE; 809 ierr = ReCompute(ts,tjsch,e->stepnum,stepnum);CHKERRQ(ierr); 810 } 811 if ((stack->solution_only && e->stepnum+1 == stepnum) || (!stack->solution_only && e->stepnum == stepnum)) { 812 ierr = StackPop(stack,&e);CHKERRQ(ierr); 813 ierr = ElementDestroy(stack,e);CHKERRQ(ierr); 814 } 815 if (tjsch->rctx->reverseonestep) tjsch->rctx->reverseonestep = PETSC_FALSE; 816 PetscFunctionReturn(0); 817 } 818 819 #undef __FUNCT__ 820 #define __FUNCT__ "SetTrajRON" 821 static PetscErrorCode SetTrajRON(TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X) 822 { 823 Stack *stack = &tjsch->stack; 824 Vec *Y; 825 PetscInt i,store; 826 PetscReal timeprev; 827 StackElement e; 828 RevolveCTX *rctx = tjsch->rctx; 829 PetscErrorCode ierr; 830 831 PetscFunctionBegin; 832 if (!stack->solution_only && stepnum == 0) PetscFunctionReturn(0); 833 if (stack->solution_only && stepnum == tjsch->total_steps) PetscFunctionReturn(0); 834 ierr = ApplyRevolve(tjsch->stype,rctx,tjsch->total_steps,stepnum,stepnum,&store);CHKERRQ(ierr); 835 if (store == 1) { 836 if (rctx->check != stack->top+1) { /* overwrite some non-top checkpoint in the stack */ 837 ierr = StackFind(stack,&e,rctx->check);CHKERRQ(ierr); 838 ierr = VecCopy(X,e->X);CHKERRQ(ierr); 839 if (stack->numY > 0 && !stack->solution_only) { 840 ierr = TSGetStages(ts,&stack->numY,&Y);CHKERRQ(ierr); 841 for (i=0;i<stack->numY;i++) { 842 ierr = VecCopy(Y[i],e->Y[i]);CHKERRQ(ierr); 843 } 844 } 845 e->stepnum = stepnum; 846 e->time = time; 847 ierr = TSGetPrevTime(ts,&timeprev);CHKERRQ(ierr); 848 e->timeprev = timeprev; 849 } else { 850 if (stepnum < stack->top) SETERRQ(tjsch->comm,PETSC_ERR_MEMC,"Illegal modification of a non-top stack element"); 851 ierr = ElementCreate(ts,stack,&e,stepnum,time,X);CHKERRQ(ierr); 852 ierr = StackPush(stack,e);CHKERRQ(ierr); 853 } 854 } 855 PetscFunctionReturn(0); 856 } 857 858 #undef __FUNCT__ 859 #define __FUNCT__ "GetTrajRON" 860 static PetscErrorCode GetTrajRON(TS ts,TJScheduler *tjsch,PetscInt stepnum) 861 { 862 Stack *stack = &tjsch->stack; 863 PetscInt whattodo,shift; 864 PetscReal stepsize; 865 StackElement e; 866 PetscErrorCode ierr; 867 868 PetscFunctionBegin; 869 if (stepnum == 0 || stepnum == tjsch->total_steps) { 870 ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr); 871 ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr); 872 if (tjsch->rctx->reverseonestep) tjsch->rctx->reverseonestep = PETSC_FALSE; 873 PetscFunctionReturn(0); 874 } 875 tjsch->rctx->capo = stepnum; 876 tjsch->rctx->oldcapo = tjsch->rctx->capo; 877 shift = 0; 878 whattodo = revolve_action(&tjsch->rctx->check,&tjsch->rctx->capo,&tjsch->rctx->fine,tjsch->rctx->snaps_in,&tjsch->rctx->info,&tjsch->rctx->where); /* whattodo=restore */ 879 if (whattodo == 8) whattodo = 5; 880 printwhattodo(whattodo,tjsch->rctx,shift); 881 /* restore a checkpoint */ 882 ierr = StackFind(stack,&e,tjsch->rctx->check);CHKERRQ(ierr); 883 ierr = UpdateTS(ts,stack,e);CHKERRQ(ierr); 884 if (!stack->solution_only) { /* whattodo must be 5 */ 885 /* ask Revolve what to do next */ 886 tjsch->rctx->oldcapo = tjsch->rctx->capo; 887 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*/ 888 printwhattodo(whattodo,tjsch->rctx,shift); 889 if (whattodo == 3 || whattodo == 4) tjsch->rctx->reverseonestep = PETSC_TRUE; 890 if (whattodo == 1) tjsch->rctx->stepsleft = tjsch->rctx->capo-tjsch->rctx->oldcapo; 891 if (!tjsch->rctx->reverseonestep && tjsch->rctx->stepsleft > 0) { 892 tjsch->rctx->stepsleft--; 893 PetscPrintf(PETSC_COMM_WORLD,"\x1B[35mSkip the step from %D to %D (already checkpointed)\033[0m\n",tjsch->rctx->oldcapo,tjsch->rctx->oldcapo+1); 894 } 895 } 896 if (stack->solution_only || (!stack->solution_only && e->stepnum < stepnum)) { 897 tjsch->recompute = PETSC_TRUE; 898 ierr = ReCompute(ts,tjsch,e->stepnum,stepnum);CHKERRQ(ierr); 899 } 900 if (tjsch->rctx->reverseonestep) tjsch->rctx->reverseonestep = PETSC_FALSE; 901 PetscFunctionReturn(0); 902 } 903 904 #undef __FUNCT__ 905 #define __FUNCT__ "SetTrajTLR" 906 static PetscErrorCode SetTrajTLR(TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X) 907 { 908 Stack *stack = &tjsch->stack; 909 PetscInt store,localstepnum,laststridesize; 910 StackElement e; 911 RevolveCTX *rctx = tjsch->rctx; 912 PetscBool resetrevolve = PETSC_FALSE; 913 PetscErrorCode ierr; 914 915 PetscFunctionBegin; 916 if (!stack->solution_only && stepnum == 0) PetscFunctionReturn(0); 917 if (stack->solution_only && stepnum == tjsch->total_steps) PetscFunctionReturn(0); 918 919 localstepnum = stepnum%tjsch->stride; 920 laststridesize = tjsch->total_steps%tjsch->stride; 921 if (!laststridesize) laststridesize = tjsch->stride; 922 923 if (!tjsch->recompute) { /* skip last stride */ 924 ierr = TopLevelStore(ts,tjsch,stepnum,localstepnum,laststridesize,&resetrevolve);CHKERRQ(ierr); 925 /* different starting points for last stride between solutin_only and !solutin_only */ 926 if (!stack->solution_only && !tjsch->save_stack && stepnum <= tjsch->total_steps-laststridesize) PetscFunctionReturn(0); 927 if (stack->solution_only && !tjsch->save_stack && stepnum < tjsch->total_steps-laststridesize) PetscFunctionReturn(0); 928 } 929 930 if (tjsch->save_stack && resetrevolve) { 931 revolve_reset(); 932 revolve_create_offline(tjsch->stride,tjsch->max_cps_ram); 933 rctx = tjsch->rctx; 934 rctx->check = 0; 935 rctx->capo = 0; 936 rctx->fine = tjsch->stride; 937 tjsch->rctx->reverseonestep = PETSC_FALSE; 938 PetscFunctionReturn(0); 939 } 940 ierr = ApplyRevolve(tjsch->stype,tjsch->rctx,tjsch->total_steps,stepnum,localstepnum,&store);CHKERRQ(ierr); 941 if (store == 1) { 942 if (localstepnum < stack->top) SETERRQ(tjsch->comm,PETSC_ERR_MEMC,"Illegal modification of a non-top stack element"); 943 ierr = ElementCreate(ts,stack,&e,stepnum,time,X);CHKERRQ(ierr); 944 ierr = StackPush(stack,e);CHKERRQ(ierr); 945 } 946 947 PetscFunctionReturn(0); 948 } 949 950 #undef __FUNCT__ 951 #define __FUNCT__ "GetTrajTLR" 952 static PetscErrorCode GetTrajTLR(TS ts,TJScheduler *tjsch,PetscInt stepnum) 953 { 954 Stack *stack = &tjsch->stack; 955 PetscInt whattodo,shift; 956 PetscInt localstepnum,stridenum,laststridesize,store; 957 PetscReal stepsize; 958 StackElement e; 959 PetscErrorCode ierr; 960 961 PetscFunctionBegin; 962 localstepnum = stepnum%tjsch->stride; 963 stridenum = stepnum/tjsch->stride; 964 if (stepnum == tjsch->total_steps) { 965 ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr); 966 ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr); 967 if ( tjsch->rctx->reverseonestep) tjsch->rctx->reverseonestep = PETSC_FALSE; 968 PetscFunctionReturn(0); 969 } 970 laststridesize = tjsch->total_steps%tjsch->stride; 971 if (!laststridesize) laststridesize = tjsch->stride; 972 if (stack->solution_only) { 973 /* fill stack */ 974 if (localstepnum == 0 && stepnum <= tjsch->total_steps-laststridesize) { 975 if (tjsch->save_stack) { 976 PetscPrintf(PETSC_COMM_WORLD,"\x1B[33mLoad stack from file\033[0m\n"); 977 ierr = StackLoadAll(ts,stack,stridenum);CHKERRQ(ierr); 978 revolve_reset(); 979 revolve_create_offline(tjsch->stride,tjsch->max_cps_ram); 980 tjsch->rctx->check = 0; 981 tjsch->rctx->capo = 0; 982 tjsch->rctx->fine = tjsch->stride; 983 whattodo = 0; 984 while(whattodo!=3) { /* stupid revolve */ 985 whattodo = revolve_action(&tjsch->rctx->check,&tjsch->rctx->capo,&tjsch->rctx->fine,tjsch->rctx->snaps_in,&tjsch->rctx->info,&tjsch->rctx->where); 986 } 987 tjsch->recompute = PETSC_TRUE; 988 tjsch->skip_trajectory = PETSC_TRUE; 989 ierr = ReCompute(ts,tjsch,stridenum*tjsch->stride-1,stridenum*tjsch->stride);CHKERRQ(ierr); 990 tjsch->skip_trajectory = PETSC_FALSE; 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 tjsch->recompute = PETSC_TRUE; 1000 ierr = ReCompute(ts,tjsch,(stridenum-1)*tjsch->stride,stridenum*tjsch->stride);CHKERRQ(ierr); 1001 } 1002 PetscFunctionReturn(0); 1003 } 1004 /* restore a checkpoint */ 1005 ierr = StackTop(stack,&e);CHKERRQ(ierr); 1006 ierr = UpdateTS(ts,stack,e);CHKERRQ(ierr); 1007 /* start with restoring a checkpoint */ 1008 tjsch->rctx->capo = stepnum; 1009 tjsch->rctx->oldcapo = tjsch->rctx->capo; 1010 shift = stepnum-localstepnum; 1011 whattodo = revolve_action(&tjsch->rctx->check,&tjsch->rctx->capo,&tjsch->rctx->fine,tjsch->rctx->snaps_in,&tjsch->rctx->info,&tjsch->rctx->where); 1012 printwhattodo(whattodo,tjsch->rctx,shift); 1013 tjsch->recompute = PETSC_TRUE; 1014 ierr = ReCompute(ts,tjsch,e->stepnum,stepnum);CHKERRQ(ierr); 1015 if (e->stepnum+1 == stepnum) { 1016 ierr = StackPop(stack,&e);CHKERRQ(ierr); 1017 ierr = ElementDestroy(stack,e);CHKERRQ(ierr); 1018 } 1019 } else { 1020 /* fill stack with info */ 1021 if (localstepnum == 0 && tjsch->total_steps-stepnum >= laststridesize) { 1022 if (tjsch->save_stack) { 1023 PetscPrintf(PETSC_COMM_WORLD,"\x1B[33mLoad stack from file\033[0m\n"); 1024 ierr = StackLoadAll(ts,stack,stridenum);CHKERRQ(ierr); 1025 revolve_reset(); 1026 revolve_create_offline(tjsch->stride,tjsch->max_cps_ram); 1027 tjsch->rctx->check = 0; 1028 tjsch->rctx->capo = 0; 1029 tjsch->rctx->fine = tjsch->stride; 1030 whattodo = 0; 1031 while(whattodo!=3) { /* stupid revolve */ 1032 whattodo = revolve_action(&tjsch->rctx->check,&tjsch->rctx->capo,&tjsch->rctx->fine,tjsch->rctx->snaps_in,&tjsch->rctx->info,&tjsch->rctx->where); 1033 } 1034 } else { 1035 PetscPrintf(PETSC_COMM_WORLD,"\x1B[33mLoad a single point from file\033[0m\n"); 1036 ierr = LoadSingle(ts,stack,stridenum);CHKERRQ(ierr); 1037 revolve_reset(); 1038 revolve_create_offline(tjsch->stride,tjsch->max_cps_ram); 1039 tjsch->rctx->check = 0; 1040 tjsch->rctx->capo = 0; 1041 tjsch->rctx->fine = tjsch->stride; 1042 shift = stepnum-localstepnum; 1043 whattodo = revolve_action(&tjsch->rctx->check,&tjsch->rctx->capo,&tjsch->rctx->fine,tjsch->rctx->snaps_in,&tjsch->rctx->info,&tjsch->rctx->where); 1044 printwhattodo(whattodo,tjsch->rctx,shift); 1045 ierr = ElementCreate(ts,stack,&e,(stridenum-1)*tjsch->stride+1,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 1046 ierr = StackPush(stack,e);CHKERRQ(ierr); 1047 tjsch->recompute = PETSC_TRUE; 1048 ierr = ReCompute(ts,tjsch,e->stepnum,stridenum*tjsch->stride);CHKERRQ(ierr); 1049 if ( tjsch->rctx->reverseonestep) { /* ready for the reverse step */ 1050 ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr); 1051 ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr); 1052 tjsch->rctx->reverseonestep = PETSC_FALSE; 1053 } 1054 } 1055 PetscFunctionReturn(0); 1056 } 1057 /* restore a checkpoint */ 1058 ierr = StackTop(stack,&e);CHKERRQ(ierr); 1059 ierr = UpdateTS(ts,stack,e);CHKERRQ(ierr); 1060 /* 2 revolve actions: restore a checkpoint and then advance */ 1061 ierr = ApplyRevolve(tjsch->stype,tjsch->rctx,tjsch->total_steps,stepnum,localstepnum,&store);CHKERRQ(ierr); 1062 if (!tjsch->rctx->reverseonestep && tjsch->rctx->stepsleft > 0) { 1063 tjsch->rctx->stepsleft--; 1064 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); 1065 } 1066 if (e->stepnum < stepnum) { 1067 tjsch->recompute = PETSC_TRUE; 1068 ierr = ReCompute(ts,tjsch,e->stepnum,stepnum);CHKERRQ(ierr); 1069 } 1070 if (e->stepnum == stepnum) { 1071 ierr = StackPop(stack,&e);CHKERRQ(ierr); 1072 ierr = ElementDestroy(stack,e);CHKERRQ(ierr); 1073 } 1074 } 1075 if (tjsch->rctx->reverseonestep) tjsch->rctx->reverseonestep = PETSC_FALSE; 1076 PetscFunctionReturn(0); 1077 } 1078 1079 #undef __FUNCT__ 1080 #define __FUNCT__ "SetTrajTLTR" 1081 static PetscErrorCode SetTrajTLTR(TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X) 1082 { 1083 Stack *stack = &tjsch->stack; 1084 DiskStack *diskstack = &tjsch->diskstack; 1085 PetscInt store,store2,localstepnum,stridenum,laststridesize; 1086 StackElement e; 1087 RevolveCTX *rctx = tjsch->rctx; 1088 PetscBool resetrevolve = PETSC_FALSE; 1089 PetscErrorCode ierr; 1090 1091 PetscFunctionBegin; 1092 if (stack->solution_only && stepnum == tjsch->total_steps) PetscFunctionReturn(0); 1093 1094 localstepnum = stepnum%tjsch->stride; /* index at the bottom level (inside a stride) */ 1095 stridenum = stepnum/tjsch->stride; /* index at the top level */ 1096 laststridesize = tjsch->total_steps%tjsch->stride; 1097 if (!laststridesize) laststridesize = tjsch->stride; 1098 1099 if (localstepnum == 0) { 1100 ierr = ApplyRevolve(tjsch->stype,tjsch->rctx2,(tjsch->total_steps+tjsch->stride-1)/tjsch->stride,stepnum,stridenum,&store2);CHKERRQ(ierr); 1101 } 1102 1103 if (store2 == 1 || tjsch->rctx2->reverseonestep) { 1104 ierr = TopLevelStore(ts,tjsch,stepnum,localstepnum,laststridesize,&resetrevolve); CHKERRQ(ierr); 1105 /* different starting points for last stride between solutin_only and !solutin_only */ 1106 if (!stack->solution_only && !tjsch->save_stack && stepnum <= tjsch->total_steps-laststridesize) PetscFunctionReturn(0); 1107 if (stack->solution_only && !tjsch->save_stack && stepnum < tjsch->total_steps-laststridesize) PetscFunctionReturn(0); 1108 if (tjsch->save_stack && resetrevolve) { 1109 revolve_reset(); 1110 revolve_create_offline(tjsch->stride,tjsch->max_cps_ram); 1111 rctx = tjsch->rctx; 1112 rctx->check = 0; 1113 rctx->capo = 0; 1114 rctx->fine = tjsch->stride; 1115 tjsch->rctx->reverseonestep = PETSC_FALSE; 1116 PetscFunctionReturn(0); 1117 } 1118 ierr = ApplyRevolve(tjsch->stype,tjsch->rctx,tjsch->total_steps,stepnum,localstepnum,&store);CHKERRQ(ierr); 1119 if (store == 1) { 1120 if (localstepnum < stack->top) SETERRQ(tjsch->comm,PETSC_ERR_MEMC,"Illegal modification of a non-top stack element"); 1121 ierr = ElementCreate(ts,stack,&e,stepnum,time,X);CHKERRQ(ierr); 1122 ierr = StackPush(stack,e);CHKERRQ(ierr); 1123 } 1124 } 1125 1126 /* if (!stack->solution_only && stepnum == 0) PetscFunctionReturn(0); 1127 if (stack->solution_only && stepnum == tjsch->total_steps) PetscFunctionReturn(0); 1128 1129 if (tjsch->save_stack && resetrevolve) { 1130 revolve_reset(); 1131 revolve_create_offline(tjsch->stride,tjsch->max_cps_ram); 1132 rctx = tjsch->rctx; 1133 rctx->check = 0; 1134 rctx->capo = 0; 1135 rctx->fine = tjsch->stride; 1136 tjsch->rctx->reverseonestep = PETSC_FALSE; 1137 } 1138 ierr = ApplyRevolve(tjsch->stype,tjsch->rctx,tjsch->total_steps,stepnum,localstepnum,&store);CHKERRQ(ierr); 1139 if (store == 1) { 1140 if (localstepnum < stack->top) SETERRQ(tjsch->comm,PETSC_ERR_MEMC,"Illegal modification of a non-top stack element"); 1141 ierr = ElementCreate(ts,stack,&e,stepnum,time,X);CHKERRQ(ierr); 1142 ierr = StackPush(stack,e);CHKERRQ(ierr); 1143 }*/ 1144 PetscFunctionReturn(0); 1145 } 1146 1147 #undef __FUNCT__ 1148 #define __FUNCT__ "GetTrajTLTR" 1149 static PetscErrorCode GetTrajTLTR(TS ts,TJScheduler *tjsch,PetscInt stepnum) 1150 { 1151 Stack *stack = &tjsch->stack; 1152 PetscInt whattodo,shift; 1153 PetscInt localstepnum,stridenum,laststridesize,store; 1154 PetscReal stepsize; 1155 StackElement e; 1156 PetscErrorCode ierr; 1157 1158 PetscFunctionBegin; 1159 localstepnum = stepnum%tjsch->stride; 1160 stridenum = stepnum/tjsch->stride; 1161 if (stepnum == tjsch->total_steps) { 1162 ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr); 1163 ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr); 1164 if ( tjsch->rctx->reverseonestep) tjsch->rctx->reverseonestep = PETSC_FALSE; 1165 PetscFunctionReturn(0); 1166 } 1167 laststridesize = tjsch->total_steps%tjsch->stride; 1168 if (!laststridesize) laststridesize = tjsch->stride; 1169 if (stack->solution_only) { 1170 /* fill stack */ 1171 if (localstepnum == 0 && stepnum <= tjsch->total_steps-laststridesize) { 1172 if (tjsch->save_stack) { 1173 PetscPrintf(PETSC_COMM_WORLD,"\x1B[33mLoad stack from file\033[0m\n"); 1174 ierr = StackLoadAll(ts,stack,stridenum);CHKERRQ(ierr); 1175 revolve_reset(); 1176 revolve_create_offline(tjsch->stride,tjsch->max_cps_ram); 1177 tjsch->rctx->check = 0; 1178 tjsch->rctx->capo = 0; 1179 tjsch->rctx->fine = tjsch->stride; 1180 whattodo = 0; 1181 while(whattodo!=3) { /* stupid revolve */ 1182 whattodo = revolve_action(&tjsch->rctx->check,&tjsch->rctx->capo,&tjsch->rctx->fine,tjsch->rctx->snaps_in,&tjsch->rctx->info,&tjsch->rctx->where); 1183 } 1184 tjsch->recompute = PETSC_TRUE; 1185 tjsch->skip_trajectory = PETSC_TRUE; 1186 ierr = ReCompute(ts,tjsch,stridenum*tjsch->stride-1,stridenum*tjsch->stride);CHKERRQ(ierr); 1187 tjsch->skip_trajectory = PETSC_FALSE; 1188 } else { 1189 PetscPrintf(PETSC_COMM_WORLD,"\x1B[33mLoad a single point from file\033[0m\n"); 1190 ierr = LoadSingle(ts,stack,stridenum);CHKERRQ(ierr); 1191 revolve_reset(); 1192 revolve_create_offline(tjsch->stride,tjsch->max_cps_ram); 1193 tjsch->rctx->check = 0; 1194 tjsch->rctx->capo = 0; 1195 tjsch->rctx->fine = tjsch->stride; 1196 tjsch->recompute = PETSC_TRUE; 1197 ierr = ReCompute(ts,tjsch,(stridenum-1)*tjsch->stride,stridenum*tjsch->stride);CHKERRQ(ierr); 1198 } 1199 PetscFunctionReturn(0); 1200 } 1201 /* restore a checkpoint */ 1202 ierr = StackTop(stack,&e);CHKERRQ(ierr); 1203 ierr = UpdateTS(ts,stack,e);CHKERRQ(ierr); 1204 /* start with restoring a checkpoint */ 1205 tjsch->rctx->capo = stepnum; 1206 tjsch->rctx->oldcapo = tjsch->rctx->capo; 1207 shift = stepnum-localstepnum; 1208 whattodo = revolve_action(&tjsch->rctx->check,&tjsch->rctx->capo,&tjsch->rctx->fine,tjsch->rctx->snaps_in,&tjsch->rctx->info,&tjsch->rctx->where); 1209 printwhattodo(whattodo,tjsch->rctx,shift); 1210 tjsch->recompute = PETSC_TRUE; 1211 ierr = ReCompute(ts,tjsch,e->stepnum,stepnum);CHKERRQ(ierr); 1212 if (e->stepnum+1 == stepnum) { 1213 ierr = StackPop(stack,&e);CHKERRQ(ierr); 1214 ierr = ElementDestroy(stack,e);CHKERRQ(ierr); 1215 } 1216 } else { 1217 /* fill stack with info */ 1218 if (localstepnum == 0 && tjsch->total_steps-stepnum >= laststridesize) { 1219 if (tjsch->save_stack) { 1220 PetscPrintf(PETSC_COMM_WORLD,"\x1B[33mLoad stack from file\033[0m\n"); 1221 ierr = StackLoadAll(ts,stack,stridenum);CHKERRQ(ierr); 1222 revolve_reset(); 1223 revolve_create_offline(tjsch->stride,tjsch->max_cps_ram); 1224 tjsch->rctx->check = 0; 1225 tjsch->rctx->capo = 0; 1226 tjsch->rctx->fine = tjsch->stride; 1227 whattodo = 0; 1228 while(whattodo!=3) { /* stupid revolve */ 1229 whattodo = revolve_action(&tjsch->rctx->check,&tjsch->rctx->capo,&tjsch->rctx->fine,tjsch->rctx->snaps_in,&tjsch->rctx->info,&tjsch->rctx->where); 1230 } 1231 } else { 1232 PetscPrintf(PETSC_COMM_WORLD,"\x1B[33mLoad a single point from file\033[0m\n"); 1233 ierr = LoadSingle(ts,stack,stridenum-1);CHKERRQ(ierr); 1234 revolve_reset(); 1235 revolve_create_offline(tjsch->stride,tjsch->max_cps_ram); 1236 tjsch->rctx->check = 0; 1237 tjsch->rctx->capo = 0; 1238 tjsch->rctx->fine = tjsch->stride; 1239 shift = stepnum-localstepnum; 1240 whattodo = revolve_action(&tjsch->rctx->check,&tjsch->rctx->capo,&tjsch->rctx->fine,tjsch->rctx->snaps_in,&tjsch->rctx->info,&tjsch->rctx->where); 1241 printwhattodo(whattodo,tjsch->rctx,shift); 1242 ierr = ElementCreate(ts,stack,&e,(stridenum-1)*tjsch->stride+1,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 1243 ierr = StackPush(stack,e);CHKERRQ(ierr); 1244 tjsch->recompute = PETSC_TRUE; 1245 ierr = ReCompute(ts,tjsch,e->stepnum,stridenum*tjsch->stride);CHKERRQ(ierr); 1246 1247 if ( tjsch->rctx->reverseonestep) { /* ready for the reverse step */ 1248 ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr); 1249 ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr); 1250 tjsch->rctx->reverseonestep = PETSC_FALSE; 1251 } 1252 } 1253 PetscFunctionReturn(0); 1254 } 1255 /* restore a checkpoint */ 1256 ierr = StackTop(stack,&e);CHKERRQ(ierr); 1257 ierr = UpdateTS(ts,stack,e);CHKERRQ(ierr); 1258 /* 2 revolve actions: restore a checkpoint and then advance */ 1259 ierr = ApplyRevolve(tjsch->stype,tjsch->rctx,tjsch->total_steps,stepnum,localstepnum,&store);CHKERRQ(ierr); 1260 if (!tjsch->rctx->reverseonestep && tjsch->rctx->stepsleft > 0) { 1261 tjsch->rctx->stepsleft--; 1262 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); 1263 } 1264 if (e->stepnum < stepnum) { 1265 tjsch->recompute = PETSC_TRUE; 1266 ierr = ReCompute(ts,tjsch,e->stepnum,stepnum);CHKERRQ(ierr); 1267 } 1268 if (e->stepnum == stepnum) { 1269 ierr = StackPop(stack,&e);CHKERRQ(ierr); 1270 ierr = ElementDestroy(stack,e);CHKERRQ(ierr); 1271 } 1272 } 1273 if (tjsch->rctx->reverseonestep) tjsch->rctx->reverseonestep = PETSC_FALSE; 1274 PetscFunctionReturn(0); 1275 } 1276 1277 #undef __FUNCT__ 1278 #define __FUNCT__ "SetTrajRMS" 1279 static PetscErrorCode SetTrajRMS(TS ts,TJScheduler *tjsch,PetscInt stepnum,PetscReal time,Vec X) 1280 { 1281 Stack *stack = &tjsch->stack; 1282 PetscInt store; 1283 StackElement e; 1284 PetscErrorCode ierr; 1285 1286 PetscFunctionBegin; 1287 if (!stack->solution_only && stepnum == 0) PetscFunctionReturn(0); 1288 if (stack->solution_only && stepnum == tjsch->total_steps) PetscFunctionReturn(0); 1289 ierr = ApplyRevolve(tjsch->stype,tjsch->rctx,tjsch->total_steps,stepnum,stepnum,&store);CHKERRQ(ierr); 1290 if (store == 1){ 1291 if (stepnum < stack->top) SETERRQ(tjsch->comm,PETSC_ERR_MEMC,"Illegal modification of a non-top stack element"); 1292 ierr = ElementCreate(ts,stack,&e,stepnum,time,X);CHKERRQ(ierr); 1293 ierr = StackPush(stack,e);CHKERRQ(ierr); 1294 } else if (store == 2) { 1295 ierr = DumpSingle(ts,stack,tjsch->rctx->check+1);CHKERRQ(ierr); 1296 } 1297 PetscFunctionReturn(0); 1298 } 1299 1300 #undef __FUNCT__ 1301 #define __FUNCT__ "GetTrajRMS" 1302 static PetscErrorCode GetTrajRMS(TS ts,TJScheduler *tjsch,PetscInt stepnum) 1303 { 1304 Stack *stack = &tjsch->stack; 1305 PetscInt whattodo,shift; 1306 PetscInt restart; 1307 PetscBool ondisk; 1308 PetscReal stepsize; 1309 StackElement e; 1310 PetscErrorCode ierr; 1311 1312 PetscFunctionBegin; 1313 if (stepnum == 0 || stepnum == tjsch->total_steps) { 1314 ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr); 1315 ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr); 1316 if (tjsch->rctx->reverseonestep) tjsch->rctx->reverseonestep = PETSC_FALSE; 1317 PetscFunctionReturn(0); 1318 } 1319 tjsch->rctx->capo = stepnum; 1320 tjsch->rctx->oldcapo = tjsch->rctx->capo; 1321 shift = 0; 1322 whattodo = revolve_action(&tjsch->rctx->check,&tjsch->rctx->capo,&tjsch->rctx->fine,tjsch->rctx->snaps_in,&tjsch->rctx->info,&tjsch->rctx->where); /* whattodo=restore */ 1323 printwhattodo(whattodo,tjsch->rctx,shift); 1324 /* restore a checkpoint */ 1325 restart = tjsch->rctx->capo; 1326 if (!tjsch->rctx->where) { 1327 ondisk = PETSC_TRUE; 1328 ierr = LoadSingle(ts,stack,tjsch->rctx->check+1);CHKERRQ(ierr); 1329 ierr = TSSetTimeStep(ts,ts->ptime_prev-ts->ptime);CHKERRQ(ierr); 1330 } else { 1331 ondisk = PETSC_FALSE; 1332 ierr = StackTop(stack,&e);CHKERRQ(ierr); 1333 ierr = UpdateTS(ts,stack,e);CHKERRQ(ierr); 1334 } 1335 if (!stack->solution_only) { /* whattodo must be 5 or 8 */ 1336 /* ask Revolve what to do next */ 1337 tjsch->rctx->oldcapo = tjsch->rctx->capo; 1338 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*/ 1339 printwhattodo(whattodo,tjsch->rctx,shift); 1340 if (whattodo == 3 || whattodo == 4) tjsch->rctx->reverseonestep = PETSC_TRUE; 1341 if (whattodo == 1) tjsch->rctx->stepsleft = tjsch->rctx->capo-tjsch->rctx->oldcapo; 1342 if (!tjsch->rctx->reverseonestep && tjsch->rctx->stepsleft > 0) { 1343 tjsch->rctx->stepsleft--; 1344 PetscPrintf(PETSC_COMM_WORLD,"\x1B[35mSkip the step from %D to %D (already checkpointed)\033[0m\n",tjsch->rctx->oldcapo,tjsch->rctx->oldcapo+1); 1345 } 1346 restart++; /* skip one step */ 1347 } 1348 if (stack->solution_only || (!stack->solution_only && restart < stepnum)) { 1349 tjsch->recompute = PETSC_TRUE; 1350 ierr = ReCompute(ts,tjsch,restart,stepnum);CHKERRQ(ierr); 1351 } 1352 if (!ondisk && ( (stack->solution_only && e->stepnum+1 == stepnum) || (!stack->solution_only && e->stepnum == stepnum) )) { 1353 ierr = StackPop(stack,&e);CHKERRQ(ierr); 1354 ierr = ElementDestroy(stack,e);CHKERRQ(ierr); 1355 } 1356 if (tjsch->rctx->reverseonestep) tjsch->rctx->reverseonestep = PETSC_FALSE; 1357 PetscFunctionReturn(0); 1358 } 1359 #endif 1360 1361 #undef __FUNCT__ 1362 #define __FUNCT__ "TSTrajectorySet_Memory" 1363 PetscErrorCode TSTrajectorySet_Memory(TSTrajectory tj,TS ts,PetscInt stepnum,PetscReal time,Vec X) 1364 { 1365 TJScheduler *tjsch = (TJScheduler*)tj->data; 1366 PetscErrorCode ierr; 1367 1368 PetscFunctionBegin; 1369 if (!tjsch->recompute) { /* use global stepnum in the forward sweep */ 1370 ierr = TSGetTotalSteps(ts,&stepnum);CHKERRQ(ierr); 1371 } 1372 /* for consistency */ 1373 if (!tjsch->recompute && stepnum == 0) ts->ptime_prev = ts->ptime-ts->time_step; 1374 switch (tjsch->stype) { 1375 case NONE: 1376 ierr = SetTrajN(ts,tjsch,stepnum,time,X);CHKERRQ(ierr); 1377 break; 1378 case TWO_LEVEL_NOREVOLVE: 1379 ierr = SetTrajTLNR(ts,tjsch,stepnum,time,X);CHKERRQ(ierr); 1380 break; 1381 #ifdef PETSC_HAVE_REVOLVE 1382 case TWO_LEVEL_REVOLVE: 1383 ierr = SetTrajTLR(ts,tjsch,stepnum,time,X);CHKERRQ(ierr); 1384 break; 1385 case TWO_LEVEL_TWO_REVOLVE: 1386 ierr = SetTrajTLTR(ts,tjsch,stepnum,time,X);CHKERRQ(ierr); 1387 break; 1388 case REVOLVE_OFFLINE: 1389 ierr = SetTrajROF(ts,tjsch,stepnum,time,X);CHKERRQ(ierr); 1390 break; 1391 case REVOLVE_ONLINE: 1392 ierr = SetTrajRON(ts,tjsch,stepnum,time,X);CHKERRQ(ierr); 1393 break; 1394 case REVOLVE_MULTISTAGE: 1395 ierr = SetTrajRMS(ts,tjsch,stepnum,time,X);CHKERRQ(ierr); 1396 break; 1397 #endif 1398 default: 1399 break; 1400 } 1401 PetscFunctionReturn(0); 1402 } 1403 1404 #undef __FUNCT__ 1405 #define __FUNCT__ "TSTrajectoryGet_Memory" 1406 PetscErrorCode TSTrajectoryGet_Memory(TSTrajectory tj,TS ts,PetscInt stepnum,PetscReal *t) 1407 { 1408 TJScheduler *tjsch = (TJScheduler*)tj->data; 1409 PetscErrorCode ierr; 1410 1411 PetscFunctionBegin; 1412 ierr = TSGetTotalSteps(ts,&stepnum);CHKERRQ(ierr); 1413 if (stepnum == 0) PetscFunctionReturn(0); 1414 switch (tjsch->stype) { 1415 case NONE: 1416 ierr = GetTrajN(ts,tjsch,stepnum);CHKERRQ(ierr); 1417 break; 1418 case TWO_LEVEL_NOREVOLVE: 1419 ierr = GetTrajTLNR(ts,tjsch,stepnum);CHKERRQ(ierr); 1420 break; 1421 #ifdef PETSC_HAVE_REVOLVE 1422 case TWO_LEVEL_REVOLVE: 1423 ierr = GetTrajTLR(ts,tjsch,stepnum);CHKERRQ(ierr); 1424 break; 1425 case TWO_LEVEL_TWO_REVOLVE: 1426 ierr = GetTrajTLTR(ts,tjsch,stepnum);CHKERRQ(ierr); 1427 break; 1428 case REVOLVE_OFFLINE: 1429 ierr = GetTrajROF(ts,tjsch,stepnum);CHKERRQ(ierr); 1430 break; 1431 case REVOLVE_ONLINE: 1432 ierr = GetTrajRON(ts,tjsch,stepnum);CHKERRQ(ierr); 1433 break; 1434 case REVOLVE_MULTISTAGE: 1435 ierr = GetTrajRMS(ts,tjsch,stepnum);CHKERRQ(ierr); 1436 break; 1437 #endif 1438 default: 1439 break; 1440 } 1441 PetscFunctionReturn(0); 1442 } 1443 1444 #undef __FUNCT__ 1445 #define __FUNCT__ "TSTrajectorySetStride_Memory" 1446 PetscErrorCode TSTrajectorySetStride_Memory(TSTrajectory tj,TS ts,PetscInt stride) 1447 { 1448 TJScheduler *tjsch = (TJScheduler*)tj->data; 1449 1450 PetscFunctionBegin; 1451 tjsch->stride = stride; 1452 PetscFunctionReturn(0); 1453 } 1454 1455 #undef __FUNCT__ 1456 #define __FUNCT__ "TSTrajectorySetMaxCpsRAM_Memory" 1457 PetscErrorCode TSTrajectorySetMaxCpsRAM_Memory(TSTrajectory tj,TS ts,PetscInt max_cps_ram) 1458 { 1459 TJScheduler *tjsch = (TJScheduler*)tj->data; 1460 1461 PetscFunctionBegin; 1462 tjsch->max_cps_ram = max_cps_ram; 1463 PetscFunctionReturn(0); 1464 } 1465 1466 #undef __FUNCT__ 1467 #define __FUNCT__ "TSTrajectorySetMaxCpsDisk_Memory" 1468 PetscErrorCode TSTrajectorySetMaxCpsDisk_Memory(TSTrajectory tj,TS ts,PetscInt max_cps_disk) 1469 { 1470 TJScheduler *tjsch = (TJScheduler*)tj->data; 1471 1472 PetscFunctionBegin; 1473 tjsch->max_cps_disk = max_cps_disk; 1474 PetscFunctionReturn(0); 1475 } 1476 1477 #ifdef PETSC_HAVE_REVOLVE 1478 #undef __FUNCT__ 1479 #define __FUNCT__ "TSTrajectorySetRevolveOnline" 1480 PetscErrorCode TSTrajectorySetRevolveOnline(TSTrajectory tj,PetscBool use_online) 1481 { 1482 TJScheduler *tjsch = (TJScheduler*)tj->data; 1483 1484 PetscFunctionBegin; 1485 tjsch->use_online = use_online; 1486 PetscFunctionReturn(0); 1487 } 1488 #endif 1489 1490 #undef __FUNCT__ 1491 #define __FUNCT__ "TSTrajectorySetSaveStack" 1492 PetscErrorCode TSTrajectorySetSaveStack(TSTrajectory tj,PetscBool save_stack) 1493 { 1494 TJScheduler *tjsch = (TJScheduler*)tj->data; 1495 1496 PetscFunctionBegin; 1497 tjsch->save_stack = save_stack; 1498 PetscFunctionReturn(0); 1499 } 1500 1501 #undef __FUNCT__ 1502 #define __FUNCT__ "TSTrajectorySetSolutionOnly" 1503 PetscErrorCode TSTrajectorySetSolutionOnly(TSTrajectory tj,PetscBool solution_only) 1504 { 1505 TJScheduler *tjsch = (TJScheduler*)tj->data; 1506 Stack *stack = &tjsch->stack; 1507 1508 PetscFunctionBegin; 1509 stack->solution_only = solution_only; 1510 PetscFunctionReturn(0); 1511 } 1512 1513 #undef __FUNCT__ 1514 #define __FUNCT__ "TSTrajectorySetFromOptions_Memory" 1515 PetscErrorCode TSTrajectorySetFromOptions_Memory(PetscOptions *PetscOptionsObject,TSTrajectory tj) 1516 { 1517 TJScheduler *tjsch = (TJScheduler*)tj->data; 1518 PetscErrorCode ierr; 1519 1520 PetscFunctionBegin; 1521 ierr = PetscOptionsHead(PetscOptionsObject,"Memory based TS trajectory options");CHKERRQ(ierr); 1522 { 1523 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); 1524 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); 1525 ierr = PetscOptionsInt("-tstrajectory_stride","Stride to save checkpoints to file","TSTrajectorySetStride_Memory",tjsch->stride,&tjsch->stride,NULL);CHKERRQ(ierr); 1526 #ifdef PETSC_HAVE_REVOLVE 1527 ierr = PetscOptionsBool("-tstrajectory_revolve_online","Trick TS trajectory into using online mode of revolve","TSTrajectorySetRevolveOnline",tjsch->use_online,&tjsch->use_online,NULL);CHKERRQ(ierr); 1528 #endif 1529 ierr = PetscOptionsBool("-tstrajectory_save_stack","Save all stack to disk","TSTrajectorySetSaveStack",tjsch->save_stack,&tjsch->save_stack,NULL);CHKERRQ(ierr); 1530 ierr = PetscOptionsBool("-tstrajectory_solution_only","Checkpoint solution only","TSTrajectorySetSolutionOnly",tjsch->stack.solution_only,&tjsch->stack.solution_only,NULL);CHKERRQ(ierr); 1531 } 1532 ierr = PetscOptionsTail();CHKERRQ(ierr); 1533 PetscFunctionReturn(0); 1534 } 1535 1536 #undef __FUNCT__ 1537 #define __FUNCT__ "TSTrajectorySetUp_Memory" 1538 PetscErrorCode TSTrajectorySetUp_Memory(TSTrajectory tj,TS ts) 1539 { 1540 TJScheduler *tjsch = (TJScheduler*)tj->data; 1541 Stack *stack = &tjsch->stack; 1542 #ifdef PETSC_HAVE_REVOLVE 1543 RevolveCTX *rctx; 1544 #endif 1545 PetscInt numY,ns; 1546 PetscBool flg; 1547 PetscErrorCode ierr; 1548 1549 PetscFunctionBegin; 1550 PetscStrcmp(((PetscObject)ts->adapt)->type_name,TSADAPTNONE,&flg); 1551 if (flg) { /* fixed time step */ 1552 tjsch->total_steps = PetscMin(ts->max_steps,(PetscInt)(ceil(ts->max_time/ts->time_step))); 1553 } 1554 if (tjsch->max_cps_ram > 1) stack->stacksize = tjsch->max_cps_ram; 1555 if (tjsch->stride > 1) { /* two level mode works for both fixed time step and adaptive time step */ 1556 if (tjsch->max_cps_disk <= 1 && tjsch->max_cps_ram > 1 && tjsch->max_cps_ram < tjsch->stride-1) { /* use revolve_offline for each stride */ 1557 tjsch->stype = TWO_LEVEL_REVOLVE; 1558 } 1559 if (tjsch->max_cps_disk > 1 && tjsch->max_cps_ram > 1 && tjsch->max_cps_ram <= tjsch->stride-1) { /* use revolve_offline for each stride */ 1560 tjsch->stype = TWO_LEVEL_TWO_REVOLVE; 1561 tjsch->diskstack.stacksize = tjsch->max_cps_disk; 1562 } 1563 if (tjsch->max_cps_disk <= 1 && (tjsch->max_cps_ram >= tjsch->stride || tjsch->max_cps_ram == -1)) { /* checkpoint all for each stride */ 1564 tjsch->stype = TWO_LEVEL_NOREVOLVE; /* can also be handled by TWO_LEVEL_REVOLVE */ 1565 stack->stacksize = tjsch->stride-1; /* need tjsch->stride-1 at most */ 1566 } 1567 } else { 1568 if (flg) { /* fixed time step */ 1569 if (tjsch->max_cps_ram >= tjsch->total_steps-1 || tjsch->max_cps_ram < 1) { /* checkpoint all */ 1570 tjsch->stype = NONE; 1571 stack->stacksize = stack->solution_only ? tjsch->total_steps : tjsch->total_steps-1; 1572 } else { 1573 if (tjsch->max_cps_disk > 1) { /* disk can be used */ 1574 tjsch->stype = REVOLVE_MULTISTAGE; 1575 } else { /* memory only */ 1576 tjsch->stype = REVOLVE_OFFLINE; 1577 } 1578 } 1579 } else { /* adaptive time step */ 1580 tjsch->stype = REVOLVE_ONLINE; 1581 } 1582 #ifdef PETSC_HAVE_REVOLVE 1583 if (tjsch->use_online) { /* trick into online */ 1584 tjsch->stype = REVOLVE_ONLINE; 1585 stack->stacksize = tjsch->max_cps_ram; 1586 } 1587 #endif 1588 } 1589 1590 if (tjsch->stype > TWO_LEVEL_NOREVOLVE) { 1591 #ifndef PETSC_HAVE_REVOLVE 1592 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."); 1593 #else 1594 if (tjsch->stype == TWO_LEVEL_REVOLVE) revolve_create_offline(tjsch->stride,tjsch->max_cps_ram); 1595 else if (tjsch->stype == TWO_LEVEL_TWO_REVOLVE) { 1596 revolve_create_offline(tjsch->stride,tjsch->max_cps_ram); 1597 ns = (tjsch->total_steps+tjsch->stride-1)/tjsch->stride; 1598 revolve2_create_offline(ns,tjsch->max_cps_disk); 1599 } else if (tjsch->stype == REVOLVE_OFFLINE) revolve_create_offline(tjsch->total_steps,tjsch->max_cps_ram); 1600 else if (tjsch->stype == REVOLVE_ONLINE) revolve_create_online(tjsch->max_cps_ram); 1601 else if (tjsch->stype == REVOLVE_MULTISTAGE) revolve_create_multistage(tjsch->total_steps,tjsch->max_cps_ram+tjsch->max_cps_disk,tjsch->max_cps_ram); 1602 1603 ierr = PetscCalloc1(1,&rctx);CHKERRQ(ierr); 1604 rctx->snaps_in = tjsch->max_cps_ram; /* for theta methods snaps_in=2*max_cps_ram */ 1605 rctx->reverseonestep = PETSC_FALSE; 1606 rctx->check = 0; 1607 rctx->oldcapo = 0; 1608 rctx->capo = 0; 1609 rctx->info = 2; 1610 rctx->fine = (tjsch->stride > 1) ? tjsch->stride : tjsch->total_steps; 1611 1612 tjsch->rctx = rctx; 1613 if (tjsch->stype == REVOLVE_ONLINE) rctx->fine = -1; 1614 #endif 1615 } 1616 1617 tjsch->recompute = PETSC_FALSE; 1618 tjsch->comm = PetscObjectComm((PetscObject)ts); 1619 ierr = TSGetStages(ts,&numY,PETSC_IGNORE);CHKERRQ(ierr); 1620 ierr = StackCreate(stack,stack->stacksize,numY);CHKERRQ(ierr); 1621 PetscFunctionReturn(0); 1622 } 1623 1624 #undef __FUNCT__ 1625 #define __FUNCT__ "TSTrajectoryDestroy_Memory" 1626 PETSC_EXTERN PetscErrorCode TSTrajectoryDestroy_Memory(TSTrajectory tj) 1627 { 1628 TJScheduler *tjsch = (TJScheduler*)tj->data; 1629 PetscErrorCode ierr; 1630 1631 PetscFunctionBegin; 1632 if (tjsch->stype > TWO_LEVEL_NOREVOLVE) { 1633 #ifdef PETSC_HAVE_REVOLVE 1634 revolve_reset(); 1635 if (tjsch->stype == TWO_LEVEL_TWO_REVOLVE) revolve2_reset(); 1636 #endif 1637 } 1638 ierr = StackDestroy(&tjsch->stack);CHKERRQ(ierr); 1639 #ifdef PETSC_HAVE_REVOLVE 1640 if (tjsch->stype > TWO_LEVEL_NOREVOLVE) { 1641 ierr = PetscFree(tjsch->rctx);CHKERRQ(ierr); 1642 } 1643 #endif 1644 ierr = PetscFree(tjsch);CHKERRQ(ierr); 1645 PetscFunctionReturn(0); 1646 } 1647 1648 /*MC 1649 TSTRAJECTORYMEMORY - Stores each solution of the ODE/ADE in memory 1650 1651 Level: intermediate 1652 1653 .seealso: TSTrajectoryCreate(), TS, TSTrajectorySetType() 1654 1655 M*/ 1656 #undef __FUNCT__ 1657 #define __FUNCT__ "TSTrajectoryCreate_Memory" 1658 PETSC_EXTERN PetscErrorCode TSTrajectoryCreate_Memory(TSTrajectory tj,TS ts) 1659 { 1660 TJScheduler *tjsch; 1661 PetscErrorCode ierr; 1662 1663 PetscFunctionBegin; 1664 tj->ops->set = TSTrajectorySet_Memory; 1665 tj->ops->get = TSTrajectoryGet_Memory; 1666 tj->ops->setup = TSTrajectorySetUp_Memory; 1667 tj->ops->destroy = TSTrajectoryDestroy_Memory; 1668 tj->ops->setfromoptions = TSTrajectorySetFromOptions_Memory; 1669 1670 ierr = PetscCalloc1(1,&tjsch);CHKERRQ(ierr); 1671 tjsch->stype = NONE; 1672 tjsch->max_cps_ram = -1; /* -1 indicates that it is not set */ 1673 tjsch->max_cps_disk = -1; /* -1 indicates that it is not set */ 1674 tjsch->stride = 0; /* if not zero, two-level checkpointing will be used */ 1675 #ifdef PETSC_HAVE_REVOLVE 1676 tjsch->use_online = PETSC_FALSE; 1677 #endif 1678 tjsch->save_stack = PETSC_TRUE; 1679 1680 tjsch->stack.solution_only = PETSC_TRUE; 1681 1682 tj->data = tjsch; 1683 1684 PetscFunctionReturn(0); 1685 } 1686