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