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