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 PetscErrorCode ierr; 203 204 PetscFunctionBegin; 205 ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr); 206 ierr = TSAdaptCandidatesClear(ts->adapt);CHKERRQ(ierr); 207 PetscFunctionReturn(0); 208 } 209 210 #undef __FUNCT__ 211 #define __FUNCT__ "TSStep_SSP" 212 static PetscErrorCode TSStep_SSP(TS ts) 213 { 214 TS_SSP *ssp = (TS_SSP*)ts->data; 215 Vec sol = ts->vec_sol; 216 PetscBool stageok,accept = PETSC_TRUE; 217 PetscReal next_time_step = ts->time_step; 218 PetscErrorCode ierr; 219 220 PetscFunctionBegin; 221 ierr = (*ssp->onestep)(ts,ts->ptime,ts->time_step,sol);CHKERRQ(ierr); 222 ierr = TSPostStage(ts,ts->ptime,0,&sol);CHKERRQ(ierr); 223 ierr = TSAdaptCheckStage(ts->adapt,ts,ts->ptime+ts->time_step,sol,&stageok);CHKERRQ(ierr); 224 if(!stageok) {ts->reason = TS_DIVERGED_STEP_REJECTED; PetscFunctionReturn(0);} 225 226 ierr = TSAdaptChoose(ts->adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr); 227 if (!accept) {ts->reason = TS_DIVERGED_STEP_REJECTED; PetscFunctionReturn(0);} 228 229 ts->ptime += ts->time_step; 230 ts->time_step = next_time_step; 231 PetscFunctionReturn(0); 232 } 233 /*------------------------------------------------------------*/ 234 235 #undef __FUNCT__ 236 #define __FUNCT__ "TSReset_SSP" 237 static PetscErrorCode TSReset_SSP(TS ts) 238 { 239 TS_SSP *ssp = (TS_SSP*)ts->data; 240 PetscErrorCode ierr; 241 242 PetscFunctionBegin; 243 if (ssp->work) {ierr = VecDestroyVecs(ssp->nwork,&ssp->work);CHKERRQ(ierr);} 244 ssp->nwork = 0; 245 ssp->workout = PETSC_FALSE; 246 PetscFunctionReturn(0); 247 } 248 249 #undef __FUNCT__ 250 #define __FUNCT__ "TSDestroy_SSP" 251 static PetscErrorCode TSDestroy_SSP(TS ts) 252 { 253 TS_SSP *ssp = (TS_SSP*)ts->data; 254 PetscErrorCode ierr; 255 256 PetscFunctionBegin; 257 ierr = TSReset_SSP(ts);CHKERRQ(ierr); 258 ierr = PetscFree(ssp->type_name);CHKERRQ(ierr); 259 ierr = PetscFree(ts->data);CHKERRQ(ierr); 260 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPGetType_C",NULL);CHKERRQ(ierr); 261 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPSetType_C",NULL);CHKERRQ(ierr); 262 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPGetNumStages_C",NULL);CHKERRQ(ierr); 263 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPSetNumStages_C",NULL);CHKERRQ(ierr); 264 PetscFunctionReturn(0); 265 } 266 /*------------------------------------------------------------*/ 267 268 #undef __FUNCT__ 269 #define __FUNCT__ "TSSSPSetType" 270 /*@C 271 TSSSPSetType - set the SSP time integration scheme to use 272 273 Logically Collective 274 275 Input Arguments: 276 ts - time stepping object 277 type - type of scheme to use 278 279 Options Database Keys: 280 -ts_ssp_type <rks2>: Type of SSP method (one of) rks2 rks3 rk104 281 -ts_ssp_nstages <5>: Number of stages 282 283 Level: beginner 284 285 .seealso: TSSSP, TSSSPGetType(), TSSSPSetNumStages(), TSSSPRKS2, TSSSPRKS3, TSSSPRK104 286 @*/ 287 PetscErrorCode TSSSPSetType(TS ts,TSSSPType type) 288 { 289 PetscErrorCode ierr; 290 291 PetscFunctionBegin; 292 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 293 ierr = PetscTryMethod(ts,"TSSSPSetType_C",(TS,TSSSPType),(ts,type));CHKERRQ(ierr); 294 PetscFunctionReturn(0); 295 } 296 297 #undef __FUNCT__ 298 #define __FUNCT__ "TSSSPGetType" 299 /*@C 300 TSSSPGetType - get the SSP time integration scheme 301 302 Logically Collective 303 304 Input Argument: 305 ts - time stepping object 306 307 Output Argument: 308 type - type of scheme being used 309 310 Level: beginner 311 312 .seealso: TSSSP, TSSSPSettype(), TSSSPSetNumStages(), TSSSPRKS2, TSSSPRKS3, TSSSPRK104 313 @*/ 314 PetscErrorCode TSSSPGetType(TS ts,TSSSPType *type) 315 { 316 PetscErrorCode ierr; 317 318 PetscFunctionBegin; 319 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 320 ierr = PetscUseMethod(ts,"TSSSPGetType_C",(TS,TSSSPType*),(ts,type));CHKERRQ(ierr); 321 PetscFunctionReturn(0); 322 } 323 324 #undef __FUNCT__ 325 #define __FUNCT__ "TSSSPSetNumStages" 326 /*@ 327 TSSSPSetNumStages - set the number of stages to use with the SSP method 328 329 Logically Collective 330 331 Input Arguments: 332 ts - time stepping object 333 nstages - number of stages 334 335 Options Database Keys: 336 -ts_ssp_type <rks2>: NumStages of SSP method (one of) rks2 rks3 rk104 337 -ts_ssp_nstages <5>: Number of stages 338 339 Level: beginner 340 341 .seealso: TSSSP, TSSSPGetNumStages(), TSSSPSetNumStages(), TSSSPRKS2, TSSSPRKS3, TSSSPRK104 342 @*/ 343 PetscErrorCode TSSSPSetNumStages(TS ts,PetscInt nstages) 344 { 345 PetscErrorCode ierr; 346 347 PetscFunctionBegin; 348 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 349 ierr = PetscTryMethod(ts,"TSSSPSetNumStages_C",(TS,PetscInt),(ts,nstages));CHKERRQ(ierr); 350 PetscFunctionReturn(0); 351 } 352 353 #undef __FUNCT__ 354 #define __FUNCT__ "TSSSPGetNumStages" 355 /*@ 356 TSSSPGetNumStages - get the number of stages in the SSP time integration scheme 357 358 Logically Collective 359 360 Input Argument: 361 ts - time stepping object 362 363 Output Argument: 364 nstages - number of stages 365 366 Level: beginner 367 368 .seealso: TSSSP, TSSSPGetType(), TSSSPSetNumStages(), TSSSPRKS2, TSSSPRKS3, TSSSPRK104 369 @*/ 370 PetscErrorCode TSSSPGetNumStages(TS ts,PetscInt *nstages) 371 { 372 PetscErrorCode ierr; 373 374 PetscFunctionBegin; 375 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 376 ierr = PetscUseMethod(ts,"TSSSPGetNumStages_C",(TS,PetscInt*),(ts,nstages));CHKERRQ(ierr); 377 PetscFunctionReturn(0); 378 } 379 380 #undef __FUNCT__ 381 #define __FUNCT__ "TSSSPSetType_SSP" 382 static PetscErrorCode TSSSPSetType_SSP(TS ts,TSSSPType type) 383 { 384 PetscErrorCode ierr,(*r)(TS,PetscReal,PetscReal,Vec); 385 TS_SSP *ssp = (TS_SSP*)ts->data; 386 387 PetscFunctionBegin; 388 ierr = PetscFunctionListFind(TSSSPList,type,&r);CHKERRQ(ierr); 389 if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TS_SSP type %s given",type); 390 ssp->onestep = r; 391 ierr = PetscFree(ssp->type_name);CHKERRQ(ierr); 392 ierr = PetscStrallocpy(type,&ssp->type_name);CHKERRQ(ierr); 393 PetscFunctionReturn(0); 394 } 395 #undef __FUNCT__ 396 #define __FUNCT__ "TSSSPGetType_SSP" 397 static PetscErrorCode TSSSPGetType_SSP(TS ts,TSSSPType *type) 398 { 399 TS_SSP *ssp = (TS_SSP*)ts->data; 400 401 PetscFunctionBegin; 402 *type = ssp->type_name; 403 PetscFunctionReturn(0); 404 } 405 #undef __FUNCT__ 406 #define __FUNCT__ "TSSSPSetNumStages_SSP" 407 static PetscErrorCode TSSSPSetNumStages_SSP(TS ts,PetscInt nstages) 408 { 409 TS_SSP *ssp = (TS_SSP*)ts->data; 410 411 PetscFunctionBegin; 412 ssp->nstages = nstages; 413 PetscFunctionReturn(0); 414 } 415 #undef __FUNCT__ 416 #define __FUNCT__ "TSSSPGetNumStages_SSP" 417 static PetscErrorCode TSSSPGetNumStages_SSP(TS ts,PetscInt *nstages) 418 { 419 TS_SSP *ssp = (TS_SSP*)ts->data; 420 421 PetscFunctionBegin; 422 *nstages = ssp->nstages; 423 PetscFunctionReturn(0); 424 } 425 426 #undef __FUNCT__ 427 #define __FUNCT__ "TSSetFromOptions_SSP" 428 static PetscErrorCode TSSetFromOptions_SSP(PetscOptionItems *PetscOptionsObject,TS ts) 429 { 430 char tname[256] = TSSSPRKS2; 431 TS_SSP *ssp = (TS_SSP*)ts->data; 432 PetscErrorCode ierr; 433 PetscBool flg; 434 435 PetscFunctionBegin; 436 ierr = PetscOptionsHead(PetscOptionsObject,"SSP ODE solver options");CHKERRQ(ierr); 437 { 438 ierr = PetscOptionsFList("-ts_ssp_type","Type of SSP method","TSSSPSetType",TSSSPList,tname,tname,sizeof(tname),&flg);CHKERRQ(ierr); 439 if (flg) { 440 ierr = TSSSPSetType(ts,tname);CHKERRQ(ierr); 441 } 442 ierr = PetscOptionsInt("-ts_ssp_nstages","Number of stages","TSSSPSetNumStages",ssp->nstages,&ssp->nstages,NULL);CHKERRQ(ierr); 443 } 444 ierr = PetscOptionsTail();CHKERRQ(ierr); 445 PetscFunctionReturn(0); 446 } 447 448 #undef __FUNCT__ 449 #define __FUNCT__ "TSView_SSP" 450 static PetscErrorCode TSView_SSP(TS ts,PetscViewer viewer) 451 { 452 TS_SSP *ssp = (TS_SSP*)ts->data; 453 PetscBool ascii; 454 PetscErrorCode ierr; 455 456 PetscFunctionBegin; 457 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&ascii);CHKERRQ(ierr); 458 if (ascii) {ierr = PetscViewerASCIIPrintf(viewer," Scheme: %s\n",ssp->type_name);CHKERRQ(ierr);} 459 PetscFunctionReturn(0); 460 } 461 462 /* ------------------------------------------------------------ */ 463 464 /*MC 465 TSSSP - Explicit strong stability preserving ODE solver 466 467 Most hyperbolic conservation laws have exact solutions that are total variation diminishing (TVD) or total variation 468 bounded (TVB) although these solutions often contain discontinuities. Spatial discretizations such as Godunov's 469 scheme and high-resolution finite volume methods (TVD limiters, ENO/WENO) are designed to preserve these properties, 470 but they are usually formulated using a forward Euler time discretization or by coupling the space and time 471 discretization as in the classical Lax-Wendroff scheme. When the space and time discretization is coupled, it is very 472 difficult to produce schemes with high temporal accuracy while preserving TVD properties. An alternative is the 473 semidiscrete formulation where we choose a spatial discretization that is TVD with forward Euler and then choose a 474 time discretization that preserves the TVD property. Such integrators are called strong stability preserving (SSP). 475 476 Let c_eff be the minimum number of function evaluations required to step as far as one step of forward Euler while 477 still being SSP. Some theoretical bounds 478 479 1. There are no explicit methods with c_eff > 1. 480 481 2. There are no explicit methods beyond order 4 (for nonlinear problems) and c_eff > 0. 482 483 3. There are no implicit methods with order greater than 1 and c_eff > 2. 484 485 This integrator provides Runge-Kutta methods of order 2, 3, and 4 with maximal values of c_eff. More stages allows 486 for larger values of c_eff which improves efficiency. These implementations are low-memory and only use 2 or 3 work 487 vectors regardless of the total number of stages, so e.g. 25-stage 3rd order methods may be an excellent choice. 488 489 Methods can be chosen with -ts_ssp_type {rks2,rks3,rk104} 490 491 rks2: Second order methods with any number s>1 of stages. c_eff = (s-1)/s 492 493 rks3: Third order methods with s=n^2 stages, n>1. c_eff = (s-n)/s 494 495 rk104: A 10-stage fourth order method. c_eff = 0.6 496 497 Level: beginner 498 499 References: 500 + 1. - Ketcheson, Highly efficient strong stability preserving Runge Kutta methods with low storage implementations, SISC, 2008. 501 - 2. - Gottlieb, Ketcheson, and Shu, High order strong stability preserving time discretizations, J Scientific Computing, 2009. 502 503 .seealso: TSCreate(), TS, TSSetType() 504 505 M*/ 506 #undef __FUNCT__ 507 #define __FUNCT__ "TSCreate_SSP" 508 PETSC_EXTERN PetscErrorCode TSCreate_SSP(TS ts) 509 { 510 TS_SSP *ssp; 511 PetscErrorCode ierr; 512 513 PetscFunctionBegin; 514 ierr = TSSSPInitializePackage();CHKERRQ(ierr); 515 516 ts->ops->setup = TSSetUp_SSP; 517 ts->ops->step = TSStep_SSP; 518 ts->ops->reset = TSReset_SSP; 519 ts->ops->destroy = TSDestroy_SSP; 520 ts->ops->setfromoptions = TSSetFromOptions_SSP; 521 ts->ops->view = TSView_SSP; 522 523 ierr = PetscNewLog(ts,&ssp);CHKERRQ(ierr); 524 ts->data = (void*)ssp; 525 526 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPGetType_C",TSSSPGetType_SSP);CHKERRQ(ierr); 527 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPSetType_C",TSSSPSetType_SSP);CHKERRQ(ierr); 528 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPGetNumStages_C",TSSSPGetNumStages_SSP);CHKERRQ(ierr); 529 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSSPSetNumStages_C",TSSSPSetNumStages_SSP);CHKERRQ(ierr); 530 531 ierr = TSSSPSetType(ts,TSSSPRKS2);CHKERRQ(ierr); 532 ssp->nstages = 5; 533 PetscFunctionReturn(0); 534 } 535 536 #undef __FUNCT__ 537 #define __FUNCT__ "TSSSPInitializePackage" 538 /*@C 539 TSSSPInitializePackage - This function initializes everything in the TSSSP package. It is called 540 from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_SSP() 541 when using static libraries. 542 543 Level: developer 544 545 .keywords: TS, TSSSP, initialize, package 546 .seealso: PetscInitialize() 547 @*/ 548 PetscErrorCode TSSSPInitializePackage(void) 549 { 550 PetscErrorCode ierr; 551 552 PetscFunctionBegin; 553 if (TSSSPPackageInitialized) PetscFunctionReturn(0); 554 TSSSPPackageInitialized = PETSC_TRUE; 555 ierr = PetscFunctionListAdd(&TSSSPList,TSSSPRKS2, TSSSPStep_RK_2);CHKERRQ(ierr); 556 ierr = PetscFunctionListAdd(&TSSSPList,TSSSPRKS3, TSSSPStep_RK_3);CHKERRQ(ierr); 557 ierr = PetscFunctionListAdd(&TSSSPList,TSSSPRK104,TSSSPStep_RK_10_4);CHKERRQ(ierr); 558 ierr = PetscRegisterFinalize(TSSSPFinalizePackage);CHKERRQ(ierr); 559 PetscFunctionReturn(0); 560 } 561 562 #undef __FUNCT__ 563 #define __FUNCT__ "TSSSPFinalizePackage" 564 /*@C 565 TSSSPFinalizePackage - This function destroys everything in the TSSSP package. It is 566 called from PetscFinalize(). 567 568 Level: developer 569 570 .keywords: Petsc, destroy, package 571 .seealso: PetscFinalize() 572 @*/ 573 PetscErrorCode TSSSPFinalizePackage(void) 574 { 575 PetscErrorCode ierr; 576 577 PetscFunctionBegin; 578 TSSSPPackageInitialized = PETSC_FALSE; 579 ierr = PetscFunctionListDestroy(&TSSSPList);CHKERRQ(ierr); 580 PetscFunctionReturn(0); 581 } 582