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