1 /* 2 Code for Timestepping with explicit SSP. 3 */ 4 #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 5 6 PetscFunctionList TSSSPList = 0; 7 static PetscBool TSSSPPackageInitialized; 8 9 typedef struct { 10 PetscErrorCode (*onestep)(TS,PetscReal,PetscReal,Vec); 11 char *type_name; 12 PetscInt nstages; 13 Vec *work; 14 PetscInt nwork; 15 PetscBool workout; 16 } TS_SSP; 17 18 19 #undef __FUNCT__ 20 #define __FUNCT__ "TSSSPGetWorkVectors" 21 static PetscErrorCode TSSSPGetWorkVectors(TS ts,PetscInt n,Vec **work) 22 { 23 TS_SSP *ssp = (TS_SSP*)ts->data; 24 PetscErrorCode ierr; 25 26 PetscFunctionBegin; 27 if (ssp->workout) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Work vectors already gotten"); 28 if (ssp->nwork < n) { 29 if (ssp->nwork > 0) { 30 ierr = VecDestroyVecs(ssp->nwork,&ssp->work);CHKERRQ(ierr); 31 } 32 ierr = VecDuplicateVecs(ts->vec_sol,n,&ssp->work);CHKERRQ(ierr); 33 ssp->nwork = n; 34 } 35 *work = ssp->work; 36 ssp->workout = PETSC_TRUE; 37 PetscFunctionReturn(0); 38 } 39 40 #undef __FUNCT__ 41 #define __FUNCT__ "TSSSPRestoreWorkVectors" 42 static PetscErrorCode TSSSPRestoreWorkVectors(TS ts,PetscInt n,Vec **work) 43 { 44 TS_SSP *ssp = (TS_SSP*)ts->data; 45 46 PetscFunctionBegin; 47 if (!ssp->workout) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Work vectors have not been gotten"); 48 if (*work != ssp->work) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Wrong work vectors checked out"); 49 ssp->workout = PETSC_FALSE; 50 *work = NULL; 51 PetscFunctionReturn(0); 52 } 53 54 #undef __FUNCT__ 55 #define __FUNCT__ "TSSSPStep_RK_2" 56 /*MC 57 TSSSPRKS2 - Optimal second order SSP Runge-Kutta method, low-storage, c_eff=(s-1)/s 58 59 Pseudocode 2 of Ketcheson 2008 60 61 Level: beginner 62 63 .seealso: TSSSP, TSSSPSetType(), TSSSPSetNumStages() 64 M*/ 65 static PetscErrorCode TSSSPStep_RK_2(TS ts,PetscReal t0,PetscReal dt,Vec sol) 66 { 67 TS_SSP *ssp = (TS_SSP*)ts->data; 68 Vec *work,F; 69 PetscInt i,s; 70 PetscErrorCode ierr; 71 72 PetscFunctionBegin; 73 s = ssp->nstages; 74 ierr = TSSSPGetWorkVectors(ts,2,&work);CHKERRQ(ierr); 75 F = work[1]; 76 ierr = VecCopy(sol,work[0]);CHKERRQ(ierr); 77 for (i=0; i<s-1; i++) { 78 PetscReal stage_time = t0+dt*(i/(s-1.)); 79 ierr = TSPreStage(ts,stage_time);CHKERRQ(ierr); 80 ierr = TSComputeRHSFunction(ts,stage_time,work[0],F);CHKERRQ(ierr); 81 ierr = VecAXPY(work[0],dt/(s-1.),F);CHKERRQ(ierr); 82 } 83 ierr = TSComputeRHSFunction(ts,t0+dt,work[0],F);CHKERRQ(ierr); 84 ierr = VecAXPBYPCZ(sol,(s-1.)/s,dt/s,1./s,work[0],F);CHKERRQ(ierr); 85 ierr = TSSSPRestoreWorkVectors(ts,2,&work);CHKERRQ(ierr); 86 PetscFunctionReturn(0); 87 } 88 89 #undef __FUNCT__ 90 #define __FUNCT__ "TSSSPStep_RK_3" 91 /*MC 92 TSSSPRKS3 - Optimal third order SSP Runge-Kutta, low-storage, c_eff=(PetscSqrtReal(s)-1)/PetscSqrtReal(s), where PetscSqrtReal(s) is an integer 93 94 Pseudocode 2 of Ketcheson 2008 95 96 Level: beginner 97 98 .seealso: TSSSP, TSSSPSetType(), TSSSPSetNumStages() 99 M*/ 100 static PetscErrorCode TSSSPStep_RK_3(TS ts,PetscReal t0,PetscReal dt,Vec sol) 101 { 102 TS_SSP *ssp = (TS_SSP*)ts->data; 103 Vec *work,F; 104 PetscInt i,s,n,r; 105 PetscReal c,stage_time; 106 PetscErrorCode ierr; 107 108 PetscFunctionBegin; 109 s = ssp->nstages; 110 n = (PetscInt)(PetscSqrtReal((PetscReal)s)+0.001); 111 r = s-n; 112 if (n*n != s) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for optimal third order schemes with %d stages, must be a square number at least 4",s); 113 ierr = TSSSPGetWorkVectors(ts,3,&work);CHKERRQ(ierr); 114 F = work[2]; 115 ierr = VecCopy(sol,work[0]);CHKERRQ(ierr); 116 for (i=0; i<(n-1)*(n-2)/2; i++) { 117 c = (i<n*(n+1)/2) ? 1.*i/(s-n) : (1.*i-n)/(s-n); 118 stage_time = t0+c*dt; 119 ierr = TSPreStage(ts,stage_time);CHKERRQ(ierr); 120 ierr = TSComputeRHSFunction(ts,stage_time,work[0],F);CHKERRQ(ierr); 121 ierr = VecAXPY(work[0],dt/r,F);CHKERRQ(ierr); 122 } 123 ierr = VecCopy(work[0],work[1]);CHKERRQ(ierr); 124 for (; i<n*(n+1)/2-1; i++) { 125 c = (i<n*(n+1)/2) ? 1.*i/(s-n) : (1.*i-n)/(s-n); 126 stage_time = t0+c*dt; 127 ierr = TSPreStage(ts,stage_time);CHKERRQ(ierr); 128 ierr = TSComputeRHSFunction(ts,stage_time,work[0],F);CHKERRQ(ierr); 129 ierr = VecAXPY(work[0],dt/r,F);CHKERRQ(ierr); 130 } 131 { 132 c = (i<n*(n+1)/2) ? 1.*i/(s-n) : (1.*i-n)/(s-n); 133 stage_time = t0+c*dt; 134 ierr = TSPreStage(ts,stage_time);CHKERRQ(ierr); 135 ierr = TSComputeRHSFunction(ts,stage_time,work[0],F);CHKERRQ(ierr); 136 ierr = VecAXPBYPCZ(work[0],1.*n/(2*n-1.),(n-1.)*dt/(r*(2*n-1)),(n-1.)/(2*n-1.),work[1],F);CHKERRQ(ierr); 137 i++; 138 } 139 for (; i<s; i++) { 140 c = (i<n*(n+1)/2) ? 1.*i/(s-n) : (1.*i-n)/(s-n); 141 stage_time = t0+c*dt; 142 ierr = TSPreStage(ts,stage_time);CHKERRQ(ierr); 143 ierr = TSComputeRHSFunction(ts,stage_time,work[0],F);CHKERRQ(ierr); 144 ierr = VecAXPY(work[0],dt/r,F);CHKERRQ(ierr); 145 } 146 ierr = VecCopy(work[0],sol);CHKERRQ(ierr); 147 ierr = TSSSPRestoreWorkVectors(ts,3,&work);CHKERRQ(ierr); 148 PetscFunctionReturn(0); 149 } 150 151 #undef __FUNCT__ 152 #define __FUNCT__ "TSSSPStep_RK_10_4" 153 /*MC 154 TSSSPRKS104 - Optimal fourth order SSP Runge-Kutta, low-storage (2N), c_eff=0.6 155 156 SSPRK(10,4), Pseudocode 3 of Ketcheson 2008 157 158 Level: beginner 159 160 .seealso: TSSSP, TSSSPSetType() 161 M*/ 162 static PetscErrorCode TSSSPStep_RK_10_4(TS ts,PetscReal t0,PetscReal dt,Vec sol) 163 { 164 const PetscReal c[10] = {0, 1./6, 2./6, 3./6, 4./6, 2./6, 3./6, 4./6, 5./6, 1}; 165 Vec *work,F; 166 PetscInt i; 167 PetscReal stage_time; 168 PetscErrorCode ierr; 169 170 PetscFunctionBegin; 171 ierr = TSSSPGetWorkVectors(ts,3,&work);CHKERRQ(ierr); 172 F = work[2]; 173 ierr = VecCopy(sol,work[0]);CHKERRQ(ierr); 174 for (i=0; i<5; i++) { 175 stage_time = t0+c[i]*dt; 176 ierr = TSPreStage(ts,stage_time);CHKERRQ(ierr); 177 ierr = TSComputeRHSFunction(ts,stage_time,work[0],F);CHKERRQ(ierr); 178 ierr = VecAXPY(work[0],dt/6,F);CHKERRQ(ierr); 179 } 180 ierr = VecAXPBYPCZ(work[1],1./25,9./25,0,sol,work[0]);CHKERRQ(ierr); 181 ierr = VecAXPBY(work[0],15,-5,work[1]);CHKERRQ(ierr); 182 for (; i<9; i++) { 183 stage_time = t0+c[i]*dt; 184 ierr = TSPreStage(ts,stage_time);CHKERRQ(ierr); 185 ierr = TSComputeRHSFunction(ts,stage_time,work[0],F);CHKERRQ(ierr); 186 ierr = VecAXPY(work[0],dt/6,F);CHKERRQ(ierr); 187 } 188 stage_time = t0+dt; 189 ierr = TSPreStage(ts,stage_time);CHKERRQ(ierr); 190 ierr = TSComputeRHSFunction(ts,stage_time,work[0],F);CHKERRQ(ierr); 191 ierr = VecAXPBYPCZ(work[1],3./5,dt/10,1,work[0],F);CHKERRQ(ierr); 192 ierr = VecCopy(work[1],sol);CHKERRQ(ierr); 193 ierr = TSSSPRestoreWorkVectors(ts,3,&work);CHKERRQ(ierr); 194 PetscFunctionReturn(0); 195 } 196 197 198 #undef __FUNCT__ 199 #define __FUNCT__ "TSSetUp_SSP" 200 static PetscErrorCode TSSetUp_SSP(TS ts) 201 { 202 203 PetscFunctionBegin; 204 PetscFunctionReturn(0); 205 } 206 207 #undef __FUNCT__ 208 #define __FUNCT__ "TSStep_SSP" 209 static PetscErrorCode TSStep_SSP(TS ts) 210 { 211 TS_SSP *ssp = (TS_SSP*)ts->data; 212 Vec sol = ts->vec_sol; 213 PetscErrorCode ierr; 214 215 PetscFunctionBegin; 216 ierr = (*ssp->onestep)(ts,ts->ptime,ts->time_step,sol);CHKERRQ(ierr); 217 ts->ptime += ts->time_step; 218 ts->steps++; 219 PetscFunctionReturn(0); 220 } 221 /*------------------------------------------------------------*/ 222 #undef __FUNCT__ 223 #define __FUNCT__ "TSReset_SSP" 224 static PetscErrorCode TSReset_SSP(TS ts) 225 { 226 TS_SSP *ssp = (TS_SSP*)ts->data; 227 PetscErrorCode ierr; 228 229 PetscFunctionBegin; 230 if (ssp->work) {ierr = VecDestroyVecs(ssp->nwork,&ssp->work);CHKERRQ(ierr);} 231 ssp->nwork = 0; 232 ssp->workout = PETSC_FALSE; 233 PetscFunctionReturn(0); 234 } 235 236 #undef __FUNCT__ 237 #define __FUNCT__ "TSDestroy_SSP" 238 static PetscErrorCode TSDestroy_SSP(TS ts) 239 { 240 TS_SSP *ssp = (TS_SSP*)ts->data; 241 PetscErrorCode ierr; 242 243 PetscFunctionBegin; 244 ierr = TSReset_SSP(ts);CHKERRQ(ierr); 245 ierr = PetscFree(ssp->type_name);CHKERRQ(ierr); 246 ierr = PetscFree(ts->data);CHKERRQ(ierr); 247 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPGetType_C",NULL);CHKERRQ(ierr); 248 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPSetType_C",NULL);CHKERRQ(ierr); 249 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPGetNumStages_C",NULL);CHKERRQ(ierr); 250 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPSetNumStages_C",NULL);CHKERRQ(ierr); 251 PetscFunctionReturn(0); 252 } 253 /*------------------------------------------------------------*/ 254 255 #undef __FUNCT__ 256 #define __FUNCT__ "TSSSPSetType" 257 /*@C 258 TSSSPSetType - set the SSP time integration scheme to use 259 260 Logically Collective 261 262 Input Arguments: 263 ts - time stepping object 264 type - type of scheme to use 265 266 Options Database Keys: 267 -ts_ssp_type <rks2>: Type of SSP method (one of) rks2 rks3 rk104 268 -ts_ssp_nstages <5>: Number of stages 269 270 Level: beginner 271 272 .seealso: TSSSP, TSSSPGetType(), TSSSPSetNumStages(), TSSSPRKS2, TSSSPRKS3, TSSSPRK104 273 @*/ 274 PetscErrorCode TSSSPSetType(TS ts,TSSSPType type) 275 { 276 PetscErrorCode ierr; 277 278 PetscFunctionBegin; 279 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 280 ierr = PetscTryMethod(ts,"TSSSPSetType_C",(TS,TSSSPType),(ts,type));CHKERRQ(ierr); 281 PetscFunctionReturn(0); 282 } 283 284 #undef __FUNCT__ 285 #define __FUNCT__ "TSSSPGetType" 286 /*@C 287 TSSSPGetType - get the SSP time integration scheme 288 289 Logically Collective 290 291 Input Argument: 292 ts - time stepping object 293 294 Output Argument: 295 type - type of scheme being used 296 297 Level: beginner 298 299 .seealso: TSSSP, TSSSPSettype(), TSSSPSetNumStages(), TSSSPRKS2, TSSSPRKS3, TSSSPRK104 300 @*/ 301 PetscErrorCode TSSSPGetType(TS ts,TSSSPType *type) 302 { 303 PetscErrorCode ierr; 304 305 PetscFunctionBegin; 306 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 307 ierr = PetscTryMethod(ts,"TSSSPGetType_C",(TS,TSSSPType*),(ts,type));CHKERRQ(ierr); 308 PetscFunctionReturn(0); 309 } 310 311 #undef __FUNCT__ 312 #define __FUNCT__ "TSSSPSetNumStages" 313 /*@ 314 TSSSPSetNumStages - set the number of stages to use with the SSP method 315 316 Logically Collective 317 318 Input Arguments: 319 ts - time stepping object 320 nstages - number of stages 321 322 Options Database Keys: 323 -ts_ssp_type <rks2>: NumStages of SSP method (one of) rks2 rks3 rk104 324 -ts_ssp_nstages <5>: Number of stages 325 326 Level: beginner 327 328 .seealso: TSSSP, TSSSPGetNumStages(), TSSSPSetNumStages(), TSSSPRKS2, TSSSPRKS3, TSSSPRK104 329 @*/ 330 PetscErrorCode TSSSPSetNumStages(TS ts,PetscInt nstages) 331 { 332 PetscErrorCode ierr; 333 334 PetscFunctionBegin; 335 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 336 ierr = PetscTryMethod(ts,"TSSSPSetNumStages_C",(TS,PetscInt),(ts,nstages));CHKERRQ(ierr); 337 PetscFunctionReturn(0); 338 } 339 340 #undef __FUNCT__ 341 #define __FUNCT__ "TSSSPGetNumStages" 342 /*@ 343 TSSSPGetNumStages - get the number of stages in the SSP time integration scheme 344 345 Logically Collective 346 347 Input Argument: 348 ts - time stepping object 349 350 Output Argument: 351 nstages - number of stages 352 353 Level: beginner 354 355 .seealso: TSSSP, TSSSPGetType(), TSSSPSetNumStages(), TSSSPRKS2, TSSSPRKS3, TSSSPRK104 356 @*/ 357 PetscErrorCode TSSSPGetNumStages(TS ts,PetscInt *nstages) 358 { 359 PetscErrorCode ierr; 360 361 PetscFunctionBegin; 362 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 363 ierr = PetscTryMethod(ts,"TSSSPGetNumStages_C",(TS,PetscInt*),(ts,nstages));CHKERRQ(ierr); 364 PetscFunctionReturn(0); 365 } 366 367 #undef __FUNCT__ 368 #define __FUNCT__ "TSSSPSetType_SSP" 369 static PetscErrorCode TSSSPSetType_SSP(TS ts,TSSSPType type) 370 { 371 PetscErrorCode ierr,(*r)(TS,PetscReal,PetscReal,Vec); 372 TS_SSP *ssp = (TS_SSP*)ts->data; 373 374 PetscFunctionBegin; 375 ierr = PetscFunctionListFind(TSSSPList,type,&r);CHKERRQ(ierr); 376 if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TS_SSP type %s given",type); 377 ssp->onestep = r; 378 ierr = PetscFree(ssp->type_name);CHKERRQ(ierr); 379 ierr = PetscStrallocpy(type,&ssp->type_name);CHKERRQ(ierr); 380 PetscFunctionReturn(0); 381 } 382 #undef __FUNCT__ 383 #define __FUNCT__ "TSSSPGetType_SSP" 384 static PetscErrorCode TSSSPGetType_SSP(TS ts,TSSSPType *type) 385 { 386 TS_SSP *ssp = (TS_SSP*)ts->data; 387 388 PetscFunctionBegin; 389 *type = ssp->type_name; 390 PetscFunctionReturn(0); 391 } 392 #undef __FUNCT__ 393 #define __FUNCT__ "TSSSPSetNumStages_SSP" 394 static PetscErrorCode TSSSPSetNumStages_SSP(TS ts,PetscInt nstages) 395 { 396 TS_SSP *ssp = (TS_SSP*)ts->data; 397 398 PetscFunctionBegin; 399 ssp->nstages = nstages; 400 PetscFunctionReturn(0); 401 } 402 #undef __FUNCT__ 403 #define __FUNCT__ "TSSSPGetNumStages_SSP" 404 static PetscErrorCode TSSSPGetNumStages_SSP(TS ts,PetscInt *nstages) 405 { 406 TS_SSP *ssp = (TS_SSP*)ts->data; 407 408 PetscFunctionBegin; 409 *nstages = ssp->nstages; 410 PetscFunctionReturn(0); 411 } 412 413 #undef __FUNCT__ 414 #define __FUNCT__ "TSSetFromOptions_SSP" 415 static PetscErrorCode TSSetFromOptions_SSP(PetscOptionItems *PetscOptionsObject,TS ts) 416 { 417 char tname[256] = TSSSPRKS2; 418 TS_SSP *ssp = (TS_SSP*)ts->data; 419 PetscErrorCode ierr; 420 PetscBool flg; 421 422 PetscFunctionBegin; 423 ierr = PetscOptionsHead(PetscOptionsObject,"SSP ODE solver options");CHKERRQ(ierr); 424 { 425 ierr = PetscOptionsFList("-ts_ssp_type","Type of SSP method","TSSSPSetType",TSSSPList,tname,tname,sizeof(tname),&flg);CHKERRQ(ierr); 426 if (flg) { 427 ierr = TSSSPSetType(ts,tname);CHKERRQ(ierr); 428 } 429 ierr = PetscOptionsInt("-ts_ssp_nstages","Number of stages","TSSSPSetNumStages",ssp->nstages,&ssp->nstages,NULL);CHKERRQ(ierr); 430 } 431 ierr = PetscOptionsTail();CHKERRQ(ierr); 432 PetscFunctionReturn(0); 433 } 434 435 #undef __FUNCT__ 436 #define __FUNCT__ "TSView_SSP" 437 static PetscErrorCode TSView_SSP(TS ts,PetscViewer viewer) 438 { 439 PetscFunctionBegin; 440 PetscFunctionReturn(0); 441 } 442 443 /* ------------------------------------------------------------ */ 444 445 /*MC 446 TSSSP - Explicit strong stability preserving ODE solver 447 448 Most hyperbolic conservation laws have exact solutions that are total variation diminishing (TVD) or total variation 449 bounded (TVB) although these solutions often contain discontinuities. Spatial discretizations such as Godunov's 450 scheme and high-resolution finite volume methods (TVD limiters, ENO/WENO) are designed to preserve these properties, 451 but they are usually formulated using a forward Euler time discretization or by coupling the space and time 452 discretization as in the classical Lax-Wendroff scheme. When the space and time discretization is coupled, it is very 453 difficult to produce schemes with high temporal accuracy while preserving TVD properties. An alternative is the 454 semidiscrete formulation where we choose a spatial discretization that is TVD with forward Euler and then choose a 455 time discretization that preserves the TVD property. Such integrators are called strong stability preserving (SSP). 456 457 Let c_eff be the minimum number of function evaluations required to step as far as one step of forward Euler while 458 still being SSP. Some theoretical bounds 459 460 1. There are no explicit methods with c_eff > 1. 461 462 2. There are no explicit methods beyond order 4 (for nonlinear problems) and c_eff > 0. 463 464 3. There are no implicit methods with order greater than 1 and c_eff > 2. 465 466 This integrator provides Runge-Kutta methods of order 2, 3, and 4 with maximal values of c_eff. More stages allows 467 for larger values of c_eff which improves efficiency. These implementations are low-memory and only use 2 or 3 work 468 vectors regardless of the total number of stages, so e.g. 25-stage 3rd order methods may be an excellent choice. 469 470 Methods can be chosen with -ts_ssp_type {rks2,rks3,rk104} 471 472 rks2: Second order methods with any number s>1 of stages. c_eff = (s-1)/s 473 474 rks3: Third order methods with s=n^2 stages, n>1. c_eff = (s-n)/s 475 476 rk104: A 10-stage fourth order method. c_eff = 0.6 477 478 Level: beginner 479 480 References: 481 + 1. - Ketcheson, Highly efficient strong stability preserving Runge Kutta methods with low storage implementations, SISC, 2008. 482 - 2. - Gottlieb, Ketcheson, and Shu, High order strong stability preserving time discretizations, J Scientific Computing, 2009. 483 484 .seealso: TSCreate(), TS, TSSetType() 485 486 M*/ 487 #undef __FUNCT__ 488 #define __FUNCT__ "TSCreate_SSP" 489 PETSC_EXTERN PetscErrorCode TSCreate_SSP(TS ts) 490 { 491 TS_SSP *ssp; 492 PetscErrorCode ierr; 493 494 PetscFunctionBegin; 495 ierr = TSSSPInitializePackage();CHKERRQ(ierr); 496 497 ts->ops->setup = TSSetUp_SSP; 498 ts->ops->step = TSStep_SSP; 499 ts->ops->reset = TSReset_SSP; 500 ts->ops->destroy = TSDestroy_SSP; 501 ts->ops->setfromoptions = TSSetFromOptions_SSP; 502 ts->ops->view = TSView_SSP; 503 504 ierr = PetscNewLog(ts,&ssp);CHKERRQ(ierr); 505 ts->data = (void*)ssp; 506 507 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPGetType_C",TSSSPGetType_SSP);CHKERRQ(ierr); 508 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPSetType_C",TSSSPSetType_SSP);CHKERRQ(ierr); 509 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPGetNumStages_C",TSSSPGetNumStages_SSP);CHKERRQ(ierr); 510 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPSetNumStages_C",TSSSPSetNumStages_SSP);CHKERRQ(ierr); 511 512 ierr = TSSSPSetType(ts,TSSSPRKS2);CHKERRQ(ierr); 513 ssp->nstages = 5; 514 PetscFunctionReturn(0); 515 } 516 517 #undef __FUNCT__ 518 #define __FUNCT__ "TSSSPInitializePackage" 519 /*@C 520 TSSSPInitializePackage - This function initializes everything in the TSSSP package. It is called 521 from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_SSP() 522 when using static libraries. 523 524 Level: developer 525 526 .keywords: TS, TSSSP, initialize, package 527 .seealso: PetscInitialize() 528 @*/ 529 PetscErrorCode TSSSPInitializePackage(void) 530 { 531 PetscErrorCode ierr; 532 533 PetscFunctionBegin; 534 if (TSSSPPackageInitialized) PetscFunctionReturn(0); 535 TSSSPPackageInitialized = PETSC_TRUE; 536 ierr = PetscFunctionListAdd(&TSSSPList,TSSSPRKS2, TSSSPStep_RK_2);CHKERRQ(ierr); 537 ierr = PetscFunctionListAdd(&TSSSPList,TSSSPRKS3, TSSSPStep_RK_3);CHKERRQ(ierr); 538 ierr = PetscFunctionListAdd(&TSSSPList,TSSSPRK104,TSSSPStep_RK_10_4);CHKERRQ(ierr); 539 ierr = PetscRegisterFinalize(TSSSPFinalizePackage);CHKERRQ(ierr); 540 PetscFunctionReturn(0); 541 } 542 543 #undef __FUNCT__ 544 #define __FUNCT__ "TSSSPFinalizePackage" 545 /*@C 546 TSSSPFinalizePackage - This function destroys everything in the TSSSP package. It is 547 called from PetscFinalize(). 548 549 Level: developer 550 551 .keywords: Petsc, destroy, package 552 .seealso: PetscFinalize() 553 @*/ 554 PetscErrorCode TSSSPFinalizePackage(void) 555 { 556 PetscErrorCode ierr; 557 558 PetscFunctionBegin; 559 TSSSPPackageInitialized = PETSC_FALSE; 560 ierr = PetscFunctionListDestroy(&TSSSPList);CHKERRQ(ierr); 561 PetscFunctionReturn(0); 562 } 563