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