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(TS ts,Stack *s,PetscInt stepnum,PetscInt localstepnum,PetscInt *whattodo) 505 { 506 #ifdef PETSC_HAVE_REVOLVE 507 PetscInt shift; 508 #endif 509 RevolveCTX *rctx; 510 PetscErrorCode ierr; 511 512 PetscFunctionBegin; 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 ierr = DumpSingle(ts,s,rctx->check);CHKERRQ(ierr); 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) { 556 PetscFunctionReturn(0); 557 } else { /* store a checkpoint to RAM and ask Revolve how many time steps to advance next */ 558 rctx->oldcapo = rctx->capo; 559 *whattodo = revolve_action(&rctx->check,&rctx->capo,&rctx->fine,rctx->snaps_in,&rctx->info,&rctx->where); /* must return 1*/ 560 printwhattodo(*whattodo,rctx,shift); 561 if (s->stype == REVOLVE_ONLINE && rctx->capo >= s->total_steps-1) { 562 revolve_turn(s->total_steps,&rctx->capo,&rctx->fine); 563 printwhattodo(*whattodo,rctx,shift); 564 } 565 rctx->stepsleft = rctx->capo-rctx->oldcapo-1; 566 } 567 #endif 568 PetscFunctionReturn(0); 569 } 570 571 #undef __FUNCT__ 572 #define __FUNCT__ "ElementCreate" 573 static PetscErrorCode ElementCreate(TS ts,Stack *s,StackElement *e,PetscInt stepnum,PetscReal time,Vec X) 574 { 575 Vec *Y; 576 PetscInt i; 577 PetscReal timeprev; 578 PetscErrorCode ierr; 579 580 PetscFunctionBegin; 581 ierr = PetscCalloc1(1,e);CHKERRQ(ierr); 582 ierr = VecDuplicate(X,&(*e)->X);CHKERRQ(ierr); 583 ierr = VecCopy(X,(*e)->X);CHKERRQ(ierr); 584 if (s->numY > 0 && !s->solution_only) { 585 ierr = TSGetStages(ts,&s->numY,&Y);CHKERRQ(ierr); 586 ierr = VecDuplicateVecs(Y[0],s->numY,&(*e)->Y);CHKERRQ(ierr); 587 for (i=0;i<s->numY;i++) { 588 ierr = VecCopy(Y[i],(*e)->Y[i]);CHKERRQ(ierr); 589 } 590 } 591 (*e)->stepnum = stepnum; 592 (*e)->time = time; 593 /* for consistency */ 594 if (stepnum == 0) { 595 (*e)->timeprev = (*e)->time - ts->time_step; 596 } else { 597 ierr = TSGetPrevTime(ts,&timeprev);CHKERRQ(ierr); 598 (*e)->timeprev = timeprev; 599 } 600 PetscFunctionReturn(0); 601 } 602 603 #undef __FUNCT__ 604 #define __FUNCT__ "ElementDestroy" 605 static PetscErrorCode ElementDestroy(Stack *s,StackElement e) 606 { 607 PetscErrorCode ierr; 608 609 PetscFunctionBegin; 610 ierr = VecDestroy(&e->X);CHKERRQ(ierr); 611 if (!s->solution_only) { 612 ierr = VecDestroyVecs(s->numY,&e->Y);CHKERRQ(ierr); 613 } 614 ierr = PetscFree(e);CHKERRQ(ierr); 615 PetscFunctionReturn(0); 616 } 617 618 #undef __FUNCT__ 619 #define __FUNCT__ "UpdateTS" 620 static PetscErrorCode UpdateTS(TS ts,Stack *s,StackElement e) 621 { 622 Vec *Y; 623 PetscInt i; 624 PetscErrorCode ierr; 625 626 PetscFunctionBegin; 627 ierr = VecCopy(e->X,ts->vec_sol);CHKERRQ(ierr); 628 if (!s->solution_only) { 629 ierr = TSGetStages(ts,&s->numY,&Y);CHKERRQ(ierr); 630 for (i=0;i<s->numY;i++) { 631 ierr = VecCopy(e->Y[i],Y[i]);CHKERRQ(ierr); 632 } 633 } 634 ierr = TSSetTimeStep(ts,e->timeprev-e->time);CHKERRQ(ierr); /* stepsize will be negative */ 635 PetscFunctionReturn(0); 636 } 637 638 #undef __FUNCT__ 639 #define __FUNCT__ "ReCompute" 640 static PetscErrorCode ReCompute(TS ts,Stack *s,StackElement e,PetscInt stepnum) 641 { 642 PetscInt i; 643 PetscErrorCode ierr; 644 645 PetscFunctionBegin; 646 /* reset ts context */ 647 ts->steps = e->stepnum; /* global stepnum */ 648 ts->ptime = e->time; 649 ts->ptime_prev = e->timeprev; 650 for (i=ts->steps;i<stepnum;i++) { /* assume fixed step size */ 651 ierr = TSTrajectorySet(ts->trajectory,ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 652 ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 653 ierr = TSStep(ts);CHKERRQ(ierr); 654 if (ts->event) { 655 ierr = TSEventMonitor(ts);CHKERRQ(ierr); 656 } 657 if (!ts->steprollback) { 658 ierr = TSPostStep(ts);CHKERRQ(ierr); 659 } 660 } 661 PetscFunctionReturn(0); 662 } 663 664 #undef __FUNCT__ 665 #define __FUNCT__ "SetTrajN" 666 static PetscErrorCode SetTrajN(TS ts,Stack *s,PetscInt stepnum,PetscReal time,Vec X) 667 { 668 StackElement e; 669 PetscErrorCode ierr; 670 671 PetscFunctionBegin; 672 if (s->solution_only) { 673 /* skip the last two steps of each stride or the whole interval */ 674 if (stepnum >= s->total_steps-1 || s->recompute) PetscFunctionReturn(0); //? 675 } else { 676 /* skip the first and the last steps of each stride or the whole interval */ 677 if (stepnum == 0 || stepnum == s->total_steps) PetscFunctionReturn(0); 678 } 679 if (stepnum < s->top) SETERRQ(s->comm,PETSC_ERR_MEMC,"Illegal modification of a non-top stack element"); 680 ierr = ElementCreate(ts,s,&e,stepnum,time,X);CHKERRQ(ierr); 681 ierr = StackPush(s,e);CHKERRQ(ierr); 682 PetscFunctionReturn(0); 683 } 684 685 #undef __FUNCT__ 686 #define __FUNCT__ "GetTrajN" 687 static PetscErrorCode GetTrajN(TS ts,Stack *s,PetscInt stepnum) 688 { 689 PetscInt steps; 690 PetscReal stepsize; 691 StackElement e; 692 PetscErrorCode ierr; 693 694 PetscFunctionBegin; 695 if (stepnum == s->total_steps) { 696 ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr); 697 ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr); 698 PetscFunctionReturn(0); 699 } 700 /* restore a checkpoint */ 701 ierr = StackTop(s,&e);CHKERRQ(ierr); 702 ierr = UpdateTS(ts,s,e);CHKERRQ(ierr); 703 if (s->solution_only) {/* recompute one step */ 704 steps = ts->steps; 705 s->recompute = PETSC_TRUE; 706 ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr); 707 ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr); 708 ierr = ReCompute(ts,s,e,stepnum);CHKERRQ(ierr); 709 ts->steps = steps; 710 ts->total_steps = stepnum; 711 ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr); 712 ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr); 713 } 714 ierr = StackPop(s,&e);CHKERRQ(ierr); 715 ierr = ElementDestroy(s,e);CHKERRQ(ierr); 716 PetscFunctionReturn(0); 717 } 718 719 #undef __FUNCT__ 720 #define __FUNCT__ "SetTrajROF" 721 static PetscErrorCode SetTrajROF(TS ts,Stack *s,PetscInt stepnum,PetscReal time,Vec X) 722 { 723 PetscInt whattodo; 724 StackElement e; 725 PetscErrorCode ierr; 726 727 PetscFunctionBegin; 728 ierr = ApplyRevolve(ts,s,stepnum,stepnum,&whattodo);CHKERRQ(ierr); 729 if (whattodo == 2){ 730 if (stepnum < s->top) SETERRQ(s->comm,PETSC_ERR_MEMC,"Illegal modification of a non-top stack element"); 731 ierr = ElementCreate(ts,s,&e,stepnum,time,X);CHKERRQ(ierr); 732 ierr = StackPush(s,e);CHKERRQ(ierr); 733 } 734 PetscFunctionReturn(0); 735 } 736 737 #undef __FUNCT__ 738 #define __FUNCT__ "GetTrajROF" 739 static PetscErrorCode GetTrajROF(TS ts,Stack *s,PetscInt stepnum) 740 { 741 #ifdef PETSC_HAVE_REVOLVE 742 PetscInt whattodo,shift; 743 #endif 744 PetscInt steps; 745 PetscReal stepsize; 746 StackElement e; 747 PetscErrorCode ierr; 748 749 PetscFunctionBegin; 750 #ifdef PETSC_HAVE_REVOLVE 751 if ( s->rctx->reverseonestep) { /* ready for the reverse step */ 752 ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr); 753 ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr); 754 s->rctx->reverseonestep = PETSC_FALSE; 755 PetscFunctionReturn(0); 756 } 757 #endif 758 if ((!s->solution_only && stepnum == 0) || stepnum == s->total_steps) { 759 ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr); 760 ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr); 761 PetscFunctionReturn(0); 762 } 763 764 if (s->solution_only) { 765 #ifdef PETSC_HAVE_REVOLVE 766 s->rctx->capo = stepnum; 767 shift = 0; 768 whattodo = revolve_action(&s->rctx->check,&s->rctx->capo,&s->rctx->fine,s->rctx->snaps_in,&s->rctx->info,&s->rctx->where); 769 printwhattodo(whattodo,s->rctx,shift); 770 #endif 771 } 772 /* restore a checkpoint */ 773 ierr = StackTop(s,&e);CHKERRQ(ierr); 774 if (e && e->stepnum >= stepnum) { 775 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); 776 } 777 ierr = UpdateTS(ts,s,e);CHKERRQ(ierr); 778 if (s->solution_only || (!s->solution_only && e->stepnum<stepnum)) { /* must recompute */ 779 steps = ts->steps; 780 s->recompute = PETSC_TRUE; 781 ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr); 782 ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr); 783 ierr = ReCompute(ts,s,e,stepnum);CHKERRQ(ierr); 784 ts->steps = steps; 785 //ts->total_steps = stepnum; 786 ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr); 787 ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr); 788 #ifdef PETSC_HAVE_REVOLVE 789 s->rctx->reverseonestep = PETSC_FALSE; 790 #endif 791 } else { 792 ierr = ElementDestroy(s,e);CHKERRQ(ierr); 793 } 794 PetscFunctionReturn(0); 795 } 796 797 #undef __FUNCT__ 798 #define __FUNCT__ "SetTrajRON" 799 static PetscErrorCode SetTrajRON(TS ts,Stack *s,PetscInt stepnum,PetscReal time,Vec X) 800 { 801 PetscInt whattodo; 802 PetscReal timeprev; 803 StackElement e; 804 RevolveCTX *rctx = s->rctx; 805 PetscErrorCode ierr; 806 807 PetscFunctionBegin; 808 ierr = ApplyRevolve(ts,s,stepnum,stepnum,&whattodo);CHKERRQ(ierr); 809 if (whattodo == 2){ 810 if (rctx->check != s->top+1) { /* overwrite some non-top checkpoint in the stack*/ 811 ierr = StackFind(s,&e,rctx->check);CHKERRQ(ierr); 812 ierr = VecCopy(X,e->X);CHKERRQ(ierr); 813 e->stepnum = stepnum; 814 e->time = time; 815 ierr = TSGetPrevTime(ts,&timeprev);CHKERRQ(ierr); 816 e->timeprev = timeprev; 817 } else { 818 if (stepnum < s->top) SETERRQ(s->comm,PETSC_ERR_MEMC,"Illegal modification of a non-top stack element"); 819 ierr = ElementCreate(ts,s,&e,stepnum,time,X);CHKERRQ(ierr); 820 ierr = StackPush(s,e);CHKERRQ(ierr); 821 } 822 } 823 PetscFunctionReturn(0); 824 } 825 826 #undef __FUNCT__ 827 #define __FUNCT__ "GetTrajRON" 828 static PetscErrorCode GetTrajRON(TS ts,Stack *s,PetscInt stepnum) 829 { 830 #ifdef PETSC_HAVE_REVOLVE 831 PetscInt whattodo,shift; 832 #endif 833 PetscInt steps; 834 PetscReal stepsize; 835 StackElement e; 836 PetscErrorCode ierr; 837 838 PetscFunctionBegin; 839 #ifdef PETSC_HAVE_REVOLVE 840 if ( s->rctx->reverseonestep) { /* ready for the reverse step */ 841 ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr); 842 ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr); 843 s->rctx->reverseonestep = PETSC_FALSE; 844 PetscFunctionReturn(0); 845 } 846 #endif 847 if ((!s->solution_only && stepnum == 0) || stepnum == s->total_steps) { 848 ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr); 849 ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr); 850 PetscFunctionReturn(0); 851 } 852 853 if (s->solution_only) { 854 #ifdef PETSC_HAVE_REVOLVE 855 s->rctx->capo = stepnum; 856 shift = 0; 857 whattodo = revolve_action(&s->rctx->check,&s->rctx->capo,&s->rctx->fine,s->rctx->snaps_in,&s->rctx->info,&s->rctx->where); 858 printwhattodo(whattodo,s->rctx,shift); 859 #endif 860 } 861 /* restore a checkpoint */ 862 ierr = StackTop(s,&e);CHKERRQ(ierr); 863 if (e && e->stepnum >= stepnum) { 864 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); 865 } 866 ierr = UpdateTS(ts,s,e);CHKERRQ(ierr); 867 if (s->solution_only || (!s->solution_only && e->stepnum<stepnum)) { /* must recompute */ 868 steps = ts->steps; 869 s->recompute = PETSC_TRUE; 870 ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr); 871 ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr); 872 ierr = ReCompute(ts,s,e,stepnum);CHKERRQ(ierr); 873 ts->steps = steps; 874 //ts->total_steps = stepnum; 875 ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr); 876 ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr); 877 #ifdef PETSC_HAVE_REVOLVE 878 s->rctx->reverseonestep = PETSC_FALSE; 879 #endif 880 } 881 PetscFunctionReturn(0); 882 } 883 884 #undef __FUNCT__ 885 #define __FUNCT__ "SetTrajTLR" 886 static PetscErrorCode SetTrajTLR(TS ts,Stack *s,PetscInt stepnum,PetscReal time,Vec X) 887 { 888 PetscInt whattodo,localstepnum,id; 889 StackElement e; 890 RevolveCTX *rctx = s->rctx; 891 PetscErrorCode ierr; 892 893 PetscFunctionBegin; 894 localstepnum = stepnum%s->stride; 895 if (localstepnum == 0 && stepnum != s->total_steps && !s->recompute) { /* save to disk */ 896 id = stepnum/s->stride; 897 if (s->save_stack) { 898 if (stepnum) { /* skip step 0 */ 899 #ifdef PETSC_HAVE_REVOLVE 900 PetscPrintf(PETSC_COMM_WORLD,"\x1B[33mDump stack to file\033[0m\n"); 901 #endif 902 ierr = StackDumpAll(ts,s,id);CHKERRQ(ierr); 903 s->top = -1; /* reset top */ 904 #ifdef PETSC_HAVE_REVOLVE 905 revolve_reset(); 906 revolve_create_offline(s->stride,s->max_cps_ram); 907 rctx = s->rctx; 908 rctx->check = 0; 909 rctx->capo = 0; 910 rctx->fine = s->stride; 911 } 912 #endif 913 } else { 914 ierr = DumpSingle(ts,s,id);CHKERRQ(ierr); 915 PetscPrintf(PETSC_COMM_WORLD,"\x1B[33mDump a single point to file\033[0m\n"); 916 } 917 } 918 /* first forward sweep only checkpoints once in each stride */ 919 if (!s->recompute && !s->save_stack) PetscFunctionReturn(0); 920 921 ierr = ApplyRevolve(ts,s,stepnum,localstepnum,&whattodo);CHKERRQ(ierr); 922 if (whattodo == 2){ 923 if (localstepnum < s->top) SETERRQ(s->comm,PETSC_ERR_MEMC,"Illegal modification of a non-top stack element"); 924 ierr = ElementCreate(ts,s,&e,localstepnum,time,X);CHKERRQ(ierr); 925 ierr = StackPush(s,e);CHKERRQ(ierr); 926 } 927 PetscFunctionReturn(0); 928 } 929 930 #undef __FUNCT__ 931 #define __FUNCT__ "GetTrajTLR" 932 static PetscErrorCode GetTrajTLR(TS ts,Stack *s,PetscInt stepnum) 933 { 934 #ifdef PETSC_HAVE_REVOLVE 935 PetscInt whattodo,shift; 936 #endif 937 PetscInt i,steps,localstepnum,id; 938 PetscReal stepsize; 939 StackElement e; 940 PetscErrorCode ierr; 941 942 PetscFunctionBegin; 943 localstepnum = stepnum%s->stride; 944 if (localstepnum == 0 && stepnum != s->total_steps) { /* load from disk */ 945 #ifdef PETSC_HAVE_REVOLVE 946 PetscPrintf(PETSC_COMM_WORLD,"\x1B[33mLoad stack from file\033[0m\n"); 947 #endif 948 if (s->save_stack) { 949 id = stepnum/s->stride; 950 ierr = StackLoadAll(ts,s,id);CHKERRQ(ierr); 951 s->top = s->stacksize-1; 952 } else { 953 id = stepnum/s->stride-1; 954 ierr = LoadSingle(ts,s,id);CHKERRQ(ierr); 955 } 956 #ifdef PETSC_HAVE_REVOLVE 957 revolve_reset(); 958 revolve_create_offline(s->stride,s->max_cps_ram); 959 s->rctx->check = 0; 960 s->rctx->capo = 0; 961 s->rctx->fine = s->stride; 962 #endif 963 if (s->save_stack) { 964 #ifdef PETSC_HAVE_REVOLVE 965 whattodo = 0; 966 while(whattodo!=3) { /* stupid revolve */ 967 whattodo = revolve_action(&s->rctx->check,&s->rctx->capo,&s->rctx->fine,s->rctx->snaps_in,&s->rctx->info,&s->rctx->where); 968 } 969 #endif 970 } else { 971 /* save ts context */ 972 steps = ts->steps; 973 ts->steps = ts->total_steps; 974 s->recompute = PETSC_TRUE; 975 for (i=ts->steps;i<stepnum;i++) { /* assume fixed step size */ 976 ierr = TSTrajectorySet(ts->trajectory,ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 977 ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 978 ierr = TSStep(ts);CHKERRQ(ierr); 979 if (ts->event) { 980 ierr = TSEventMonitor(ts);CHKERRQ(ierr); 981 } 982 if (!ts->steprollback) { 983 ierr = TSPostStep(ts);CHKERRQ(ierr); 984 } 985 } 986 ts->steps = steps; 987 ts->total_steps = stepnum; 988 } 989 } 990 #ifdef PETSC_HAVE_REVOLVE 991 if ( s->rctx->reverseonestep) { /* ready for the reverse step */ 992 ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr); 993 ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr); 994 s->rctx->reverseonestep = PETSC_FALSE; 995 PetscFunctionReturn(0); 996 } 997 #endif 998 if ((!s->solution_only && stepnum == 0) || stepnum == s->total_steps) { 999 ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr); 1000 ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr); 1001 PetscFunctionReturn(0); 1002 } 1003 1004 if (s->solution_only) { 1005 #ifdef PETSC_HAVE_REVOLVE 1006 s->rctx->capo = stepnum; 1007 shift = 0; 1008 whattodo = revolve_action(&s->rctx->check,&s->rctx->capo,&s->rctx->fine,s->rctx->snaps_in,&s->rctx->info,&s->rctx->where); 1009 printwhattodo(whattodo,s->rctx,shift); 1010 #endif 1011 } 1012 /* restore a checkpoint */ 1013 ierr = StackTop(s,&e);CHKERRQ(ierr); 1014 if (e && e->stepnum >= stepnum) { 1015 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); 1016 } 1017 ierr = UpdateTS(ts,s,e);CHKERRQ(ierr); 1018 if (s->solution_only || (!s->solution_only && e->stepnum<stepnum)) { /* must recompute */ 1019 steps = ts->steps; 1020 s->recompute = PETSC_TRUE; 1021 ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr); 1022 ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr); 1023 ierr = ReCompute(ts,s,e,stepnum);CHKERRQ(ierr); 1024 ts->steps = steps; 1025 //ts->total_steps = stepnum; 1026 ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr); 1027 ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr); 1028 #ifdef PETSC_HAVE_REVOLVE 1029 s->rctx->reverseonestep = PETSC_FALSE; 1030 #endif 1031 } else { 1032 ierr = ElementDestroy(s,e);CHKERRQ(ierr); 1033 } 1034 PetscFunctionReturn(0); 1035 } 1036 1037 #undef __FUNCT__ 1038 #define __FUNCT__ "SetTrajTLNR" 1039 static PetscErrorCode SetTrajTLNR(TS ts,Stack *s,PetscInt stepnum,PetscReal time,Vec X) 1040 { 1041 PetscInt whattodo,localstepnum,id; 1042 StackElement e; 1043 PetscErrorCode ierr; 1044 1045 PetscFunctionBegin; 1046 localstepnum = stepnum%s->stride; 1047 if (localstepnum == 0 && stepnum != s->total_steps && stepnum != 0 && !s->recompute) { 1048 id = stepnum/s->stride; 1049 if (s->save_stack) { 1050 ierr = StackDumpAll(ts,s,id);CHKERRQ(ierr); 1051 s->top = -1; /* reset top */ 1052 } else { 1053 ierr = DumpSingle(ts,s,id);CHKERRQ(ierr); 1054 } 1055 } 1056 1057 ierr = ApplyRevolve(ts,s,stepnum,localstepnum,&whattodo);CHKERRQ(ierr); 1058 if (whattodo == 2){ 1059 if (localstepnum < s->top) SETERRQ(s->comm,PETSC_ERR_MEMC,"Illegal modification of a non-top stack element"); 1060 ierr = ElementCreate(ts,s,&e,localstepnum,time,X);CHKERRQ(ierr); 1061 ierr = StackPush(s,e);CHKERRQ(ierr); 1062 } 1063 PetscFunctionReturn(0); 1064 } 1065 1066 #undef __FUNCT__ 1067 #define __FUNCT__ "GetTrajTLNR" 1068 static PetscErrorCode GetTrajTLNR(TS ts,Stack *s,PetscInt stepnum) 1069 { 1070 PetscInt steps,id,localstepnum; 1071 PetscReal stepsize; 1072 StackElement e; 1073 PetscErrorCode ierr; 1074 1075 PetscFunctionBegin; 1076 localstepnum = stepnum%s->stride; 1077 if (stepnum != s->total_steps && localstepnum==0) { 1078 id = stepnum/s->stride; 1079 if (s->save_stack) { 1080 ierr = StackLoadAll(ts,s,id);CHKERRQ(ierr); 1081 } else { 1082 ierr = LoadSingle(ts,s,id);CHKERRQ(ierr); 1083 } 1084 } 1085 if ((!s->solution_only && stepnum == 0) || stepnum == s->total_steps) { 1086 ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr); 1087 ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr); 1088 PetscFunctionReturn(0); 1089 } 1090 1091 /* restore a checkpoint */ 1092 ierr = StackTop(s,&e);CHKERRQ(ierr); 1093 if (e && e->stepnum >= stepnum) { 1094 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); 1095 } 1096 ierr = UpdateTS(ts,s,e);CHKERRQ(ierr); 1097 if (s->solution_only || (!s->solution_only && e->stepnum<stepnum)) { /* must recompute */ 1098 steps = ts->steps; 1099 s->recompute = PETSC_TRUE; 1100 ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr); 1101 ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr); 1102 ierr = ReCompute(ts,s,e,stepnum);CHKERRQ(ierr); 1103 ts->steps = steps; 1104 //ts->total_steps = stepnum; 1105 ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr); 1106 ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr); 1107 } else { 1108 ierr = ElementDestroy(s,e);CHKERRQ(ierr); 1109 } 1110 PetscFunctionReturn(0); 1111 } 1112 1113 #undef __FUNCT__ 1114 #define __FUNCT__ "SetTrajRMS" 1115 static PetscErrorCode SetTrajRMS(TS ts,Stack *s,PetscInt stepnum,PetscReal time,Vec X) 1116 { 1117 PetscInt whattodo; 1118 StackElement e; 1119 PetscErrorCode ierr; 1120 1121 PetscFunctionBegin; 1122 ierr = ApplyRevolve(ts,s,stepnum,stepnum,&whattodo);CHKERRQ(ierr); 1123 if (whattodo == 2){ 1124 if (stepnum < s->top) SETERRQ(s->comm,PETSC_ERR_MEMC,"Illegal modification of a non-top stack element"); 1125 ierr = ElementCreate(ts,s,&e,stepnum,time,X);CHKERRQ(ierr); 1126 ierr = StackPush(s,e);CHKERRQ(ierr); 1127 } 1128 PetscFunctionReturn(0); 1129 } 1130 1131 #undef __FUNCT__ 1132 #define __FUNCT__ "GetTrajRMS" 1133 static PetscErrorCode GetTrajRMS(TS ts,Stack *s,PetscInt stepnum) 1134 { 1135 #ifdef PETSC_HAVE_REVOLVE 1136 PetscInt whattodo,shift; 1137 #endif 1138 PetscInt steps; 1139 PetscReal stepsize; 1140 StackElement e; 1141 PetscErrorCode ierr; 1142 1143 PetscFunctionBegin; 1144 #ifdef PETSC_HAVE_REVOLVE 1145 if ( s->rctx->reverseonestep) { /* ready for the reverse step */ 1146 ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr); 1147 ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr); 1148 s->rctx->reverseonestep = PETSC_FALSE; 1149 PetscFunctionReturn(0); 1150 } 1151 #endif 1152 if ((!s->solution_only && stepnum == 0) || stepnum == s->total_steps) { 1153 ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr); 1154 ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr); 1155 PetscFunctionReturn(0); 1156 } 1157 1158 if (s->solution_only) { 1159 #ifdef PETSC_HAVE_REVOLVE 1160 s->rctx->capo = stepnum; 1161 shift = 0; 1162 whattodo = revolve_action(&s->rctx->check,&s->rctx->capo,&s->rctx->fine,s->rctx->snaps_in,&s->rctx->info,&s->rctx->where); 1163 printwhattodo(whattodo,s->rctx,shift); 1164 #endif 1165 } 1166 /* restore a checkpoint */ 1167 if (!s->rctx->where) { 1168 ierr = LoadSingle(ts,s,stepnum);CHKERRQ(ierr); 1169 } else { 1170 ierr = StackTop(s,&e);CHKERRQ(ierr); 1171 if (e && e->stepnum >= stepnum) { 1172 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); 1173 } 1174 ierr = UpdateTS(ts,s,e);CHKERRQ(ierr); 1175 } 1176 1177 if (s->solution_only || (!s->solution_only && e->stepnum<stepnum)) { /* must recompute */ 1178 steps = ts->steps; 1179 s->recompute = PETSC_TRUE; 1180 ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr); 1181 ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr); 1182 ierr = ReCompute(ts,s,e,stepnum);CHKERRQ(ierr); 1183 ts->steps = steps; 1184 //ts->total_steps = stepnum; 1185 ierr = TSGetTimeStep(ts,&stepsize);CHKERRQ(ierr); 1186 ierr = TSSetTimeStep(ts,-stepsize);CHKERRQ(ierr); 1187 #ifdef PETSC_HAVE_REVOLVE 1188 s->rctx->reverseonestep = PETSC_FALSE; 1189 #endif 1190 } else { 1191 ierr = ElementDestroy(s,e);CHKERRQ(ierr); 1192 } 1193 PetscFunctionReturn(0); 1194 } 1195 1196 #undef __FUNCT__ 1197 #define __FUNCT__ "TSTrajectorySet_Memory" 1198 PetscErrorCode TSTrajectorySet_Memory(TSTrajectory tj,TS ts,PetscInt stepnum,PetscReal time,Vec X) 1199 { 1200 Stack *s = (Stack*)tj->data; 1201 PetscErrorCode ierr; 1202 1203 PetscFunctionBegin; 1204 if (!s->recompute) { /* use global stepnum in the forward sweep */ 1205 ierr = TSGetTotalSteps(ts,&stepnum);CHKERRQ(ierr); 1206 } 1207 /* for consistency */ 1208 if (!s->recompute && stepnum == 0) ts->ptime_prev = ts->ptime-ts->time_step; 1209 switch (s->stype) { 1210 case NONE: 1211 ierr = SetTrajN(ts,s,stepnum,time,X);CHKERRQ(ierr); 1212 break; 1213 case TWO_LEVEL_NOREVOLVE: 1214 ierr = SetTrajTLNR(ts,s,stepnum,time,X);CHKERRQ(ierr); 1215 break; 1216 case TWO_LEVEL_REVOLVE: 1217 ierr = SetTrajTLR(ts,s,stepnum,time,X);CHKERRQ(ierr); 1218 break; 1219 case REVOLVE_OFFLINE: 1220 ierr = SetTrajROF(ts,s,stepnum,time,X);CHKERRQ(ierr); 1221 break; 1222 case REVOLVE_ONLINE: 1223 ierr = SetTrajRON(ts,s,stepnum,time,X);CHKERRQ(ierr); 1224 break; 1225 case REVOLVE_MULTISTAGE: 1226 ierr = SetTrajRMS(ts,s,stepnum,time,X);CHKERRQ(ierr); 1227 break; 1228 default: 1229 break; 1230 } 1231 PetscFunctionReturn(0); 1232 } 1233 1234 #undef __FUNCT__ 1235 #define __FUNCT__ "TSTrajectoryGet_Memory" 1236 PetscErrorCode TSTrajectoryGet_Memory(TSTrajectory tj,TS ts,PetscInt stepnum,PetscReal *t) 1237 { 1238 Stack *s = (Stack*)tj->data; 1239 PetscErrorCode ierr; 1240 1241 PetscFunctionBegin; 1242 ierr = TSGetTotalSteps(ts,&stepnum);CHKERRQ(ierr); 1243 if (stepnum == 0) PetscFunctionReturn(0); 1244 switch (s->stype) { 1245 case NONE: 1246 ierr = GetTrajN(ts,s,stepnum);CHKERRQ(ierr); 1247 break; 1248 case TWO_LEVEL_NOREVOLVE: 1249 ierr = GetTrajTLNR(ts,s,stepnum);CHKERRQ(ierr); 1250 break; 1251 case TWO_LEVEL_REVOLVE: 1252 ierr = GetTrajTLR(ts,s,stepnum);CHKERRQ(ierr); 1253 break; 1254 case REVOLVE_OFFLINE: 1255 ierr = GetTrajROF(ts,s,stepnum);CHKERRQ(ierr); 1256 break; 1257 case REVOLVE_ONLINE: 1258 ierr = GetTrajRON(ts,s,stepnum);CHKERRQ(ierr); 1259 break; 1260 case REVOLVE_MULTISTAGE: 1261 ierr = GetTrajRMS(ts,s,stepnum);CHKERRQ(ierr); 1262 break; 1263 default: 1264 break; 1265 } 1266 PetscFunctionReturn(0); 1267 } 1268 1269 #undef __FUNCT__ 1270 #define __FUNCT__ "TSTrajectoryDestroy_Memory" 1271 PETSC_EXTERN PetscErrorCode TSTrajectoryDestroy_Memory(TSTrajectory tj) 1272 { 1273 Stack *s = (Stack*)tj->data; 1274 PetscErrorCode ierr; 1275 1276 PetscFunctionBegin; 1277 if (s->stype > TWO_LEVEL_NOREVOLVE) { 1278 #ifdef PETSC_HAVE_REVOLVE 1279 revolve_reset(); 1280 #endif 1281 } 1282 ierr = StackDestroy(&s);CHKERRQ(ierr); 1283 PetscFunctionReturn(0); 1284 } 1285 1286 /*MC 1287 TSTRAJECTORYMEMORY - Stores each solution of the ODE/ADE in memory 1288 1289 Level: intermediate 1290 1291 .seealso: TSTrajectoryCreate(), TS, TSTrajectorySetType() 1292 1293 M*/ 1294 #undef __FUNCT__ 1295 #define __FUNCT__ "TSTrajectoryCreate_Memory" 1296 PETSC_EXTERN PetscErrorCode TSTrajectoryCreate_Memory(TSTrajectory tj,TS ts) 1297 { 1298 Stack *s; 1299 PetscErrorCode ierr; 1300 1301 PetscFunctionBegin; 1302 tj->ops->set = TSTrajectorySet_Memory; 1303 tj->ops->get = TSTrajectoryGet_Memory; 1304 tj->ops->setup = TSTrajectorySetUp_Memory; 1305 tj->ops->destroy = TSTrajectoryDestroy_Memory; 1306 tj->ops->setfromoptions = TSTrajectorySetFromOptions_Memory; 1307 1308 ierr = PetscCalloc1(1,&s);CHKERRQ(ierr); 1309 s->stype = NONE; 1310 s->max_cps_ram = -1; /* -1 indicates that it is not set */ 1311 s->max_cps_disk = -1; /* -1 indicates that it is not set */ 1312 s->stride = 0; /* if not zero, two-level checkpointing will be used */ 1313 s->use_online = PETSC_FALSE; 1314 s->save_stack = PETSC_TRUE; 1315 s->solution_only= PETSC_TRUE; 1316 tj->data = s; 1317 PetscFunctionReturn(0); 1318 } 1319