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