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