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