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