1 /* 2 Code for timestepping with implicit generalized-\alpha method 3 for first order systems. 4 */ 5 #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 6 7 static PetscBool cited = PETSC_FALSE; 8 static const char citation[] = "@article{Jansen2000,\n" 9 " title = {A generalized-$\\alpha$ method for integrating the filtered {N}avier--{S}tokes equations with a stabilized finite element method},\n" 10 " author = {Kenneth E. Jansen and Christian H. Whiting and Gregory M. Hulbert},\n" 11 " journal = {Computer Methods in Applied Mechanics and Engineering},\n" 12 " volume = {190},\n" 13 " number = {3--4},\n" 14 " pages = {305--319},\n" 15 " year = {2000},\n" 16 " issn = {0045-7825},\n" 17 " doi = {http://dx.doi.org/10.1016/S0045-7825(00)00203-6}\n}\n"; 18 19 typedef struct { 20 PetscReal stage_time; 21 PetscReal shift_V; 22 PetscReal scale_F; 23 Vec X0, Xa, X1; 24 Vec V0, Va, V1; 25 26 PetscReal Alpha_m; 27 PetscReal Alpha_f; 28 PetscReal Gamma; 29 PetscInt order; 30 31 Vec vec_sol_prev; 32 Vec vec_lte_work; 33 34 TSStepStatus status; 35 } TS_Alpha; 36 37 static PetscErrorCode TSAlpha_StageTime(TS ts) { 38 TS_Alpha *th = (TS_Alpha *)ts->data; 39 PetscReal t = ts->ptime; 40 PetscReal dt = ts->time_step; 41 PetscReal Alpha_m = th->Alpha_m; 42 PetscReal Alpha_f = th->Alpha_f; 43 PetscReal Gamma = th->Gamma; 44 45 PetscFunctionBegin; 46 th->stage_time = t + Alpha_f * dt; 47 th->shift_V = Alpha_m / (Alpha_f * Gamma * dt); 48 th->scale_F = 1 / Alpha_f; 49 PetscFunctionReturn(0); 50 } 51 52 static PetscErrorCode TSAlpha_StageVecs(TS ts, Vec X) { 53 TS_Alpha *th = (TS_Alpha *)ts->data; 54 Vec X1 = X, V1 = th->V1; 55 Vec Xa = th->Xa, Va = th->Va; 56 Vec X0 = th->X0, V0 = th->V0; 57 PetscReal dt = ts->time_step; 58 PetscReal Alpha_m = th->Alpha_m; 59 PetscReal Alpha_f = th->Alpha_f; 60 PetscReal Gamma = th->Gamma; 61 62 PetscFunctionBegin; 63 /* V1 = 1/(Gamma*dT)*(X1-X0) + (1-1/Gamma)*V0 */ 64 PetscCall(VecWAXPY(V1, -1.0, X0, X1)); 65 PetscCall(VecAXPBY(V1, 1 - 1 / Gamma, 1 / (Gamma * dt), V0)); 66 /* Xa = X0 + Alpha_f*(X1-X0) */ 67 PetscCall(VecWAXPY(Xa, -1.0, X0, X1)); 68 PetscCall(VecAYPX(Xa, Alpha_f, X0)); 69 /* Va = V0 + Alpha_m*(V1-V0) */ 70 PetscCall(VecWAXPY(Va, -1.0, V0, V1)); 71 PetscCall(VecAYPX(Va, Alpha_m, V0)); 72 PetscFunctionReturn(0); 73 } 74 75 static PetscErrorCode TSAlpha_SNESSolve(TS ts, Vec b, Vec x) { 76 PetscInt nits, lits; 77 78 PetscFunctionBegin; 79 PetscCall(SNESSolve(ts->snes, b, x)); 80 PetscCall(SNESGetIterationNumber(ts->snes, &nits)); 81 PetscCall(SNESGetLinearSolveIterations(ts->snes, &lits)); 82 ts->snes_its += nits; 83 ts->ksp_its += lits; 84 PetscFunctionReturn(0); 85 } 86 87 /* 88 Compute a consistent initial state for the generalized-alpha method. 89 - Solve two successive backward Euler steps with halved time step. 90 - Compute the initial time derivative using backward differences. 91 - If using adaptivity, estimate the LTE of the initial step. 92 */ 93 static PetscErrorCode TSAlpha_Restart(TS ts, PetscBool *initok) { 94 TS_Alpha *th = (TS_Alpha *)ts->data; 95 PetscReal time_step; 96 PetscReal alpha_m, alpha_f, gamma; 97 Vec X0 = ts->vec_sol, X1, X2 = th->X1; 98 PetscBool stageok; 99 100 PetscFunctionBegin; 101 PetscCall(VecDuplicate(X0, &X1)); 102 103 /* Setup backward Euler with halved time step */ 104 PetscCall(TSAlphaGetParams(ts, &alpha_m, &alpha_f, &gamma)); 105 PetscCall(TSAlphaSetParams(ts, 1, 1, 1)); 106 PetscCall(TSGetTimeStep(ts, &time_step)); 107 ts->time_step = time_step / 2; 108 PetscCall(TSAlpha_StageTime(ts)); 109 th->stage_time = ts->ptime; 110 PetscCall(VecZeroEntries(th->V0)); 111 112 /* First BE step, (t0,X0) -> (t1,X1) */ 113 th->stage_time += ts->time_step; 114 PetscCall(VecCopy(X0, th->X0)); 115 PetscCall(TSPreStage(ts, th->stage_time)); 116 PetscCall(VecCopy(th->X0, X1)); 117 PetscCall(TSAlpha_SNESSolve(ts, NULL, X1)); 118 PetscCall(TSPostStage(ts, th->stage_time, 0, &X1)); 119 PetscCall(TSAdaptCheckStage(ts->adapt, ts, th->stage_time, X1, &stageok)); 120 if (!stageok) goto finally; 121 122 /* Second BE step, (t1,X1) -> (t2,X2) */ 123 th->stage_time += ts->time_step; 124 PetscCall(VecCopy(X1, th->X0)); 125 PetscCall(TSPreStage(ts, th->stage_time)); 126 PetscCall(VecCopy(th->X0, X2)); 127 PetscCall(TSAlpha_SNESSolve(ts, NULL, X2)); 128 PetscCall(TSPostStage(ts, th->stage_time, 0, &X2)); 129 PetscCall(TSAdaptCheckStage(ts->adapt, ts, th->stage_time, X2, &stageok)); 130 if (!stageok) goto finally; 131 132 /* Compute V0 ~ dX/dt at t0 with backward differences */ 133 PetscCall(VecZeroEntries(th->V0)); 134 PetscCall(VecAXPY(th->V0, -3 / ts->time_step, X0)); 135 PetscCall(VecAXPY(th->V0, +4 / ts->time_step, X1)); 136 PetscCall(VecAXPY(th->V0, -1 / ts->time_step, X2)); 137 138 /* Rough, lower-order estimate LTE of the initial step */ 139 if (th->vec_lte_work) { 140 PetscCall(VecZeroEntries(th->vec_lte_work)); 141 PetscCall(VecAXPY(th->vec_lte_work, +2, X2)); 142 PetscCall(VecAXPY(th->vec_lte_work, -4, X1)); 143 PetscCall(VecAXPY(th->vec_lte_work, +2, X0)); 144 } 145 146 finally: 147 /* Revert TSAlpha to the initial state (t0,X0) */ 148 if (initok) *initok = stageok; 149 PetscCall(TSSetTimeStep(ts, time_step)); 150 PetscCall(TSAlphaSetParams(ts, alpha_m, alpha_f, gamma)); 151 PetscCall(VecCopy(ts->vec_sol, th->X0)); 152 153 PetscCall(VecDestroy(&X1)); 154 PetscFunctionReturn(0); 155 } 156 157 static PetscErrorCode TSStep_Alpha(TS ts) { 158 TS_Alpha *th = (TS_Alpha *)ts->data; 159 PetscInt rejections = 0; 160 PetscBool stageok, accept = PETSC_TRUE; 161 PetscReal next_time_step = ts->time_step; 162 163 PetscFunctionBegin; 164 PetscCall(PetscCitationsRegister(citation, &cited)); 165 166 if (!ts->steprollback) { 167 if (th->vec_sol_prev) PetscCall(VecCopy(th->X0, th->vec_sol_prev)); 168 PetscCall(VecCopy(ts->vec_sol, th->X0)); 169 PetscCall(VecCopy(th->V1, th->V0)); 170 } 171 172 th->status = TS_STEP_INCOMPLETE; 173 while (!ts->reason && th->status != TS_STEP_COMPLETE) { 174 if (ts->steprestart) { 175 PetscCall(TSAlpha_Restart(ts, &stageok)); 176 if (!stageok) goto reject_step; 177 } 178 179 PetscCall(TSAlpha_StageTime(ts)); 180 PetscCall(VecCopy(th->X0, th->X1)); 181 PetscCall(TSPreStage(ts, th->stage_time)); 182 PetscCall(TSAlpha_SNESSolve(ts, NULL, th->X1)); 183 PetscCall(TSPostStage(ts, th->stage_time, 0, &th->Xa)); 184 PetscCall(TSAdaptCheckStage(ts->adapt, ts, th->stage_time, th->Xa, &stageok)); 185 if (!stageok) goto reject_step; 186 187 th->status = TS_STEP_PENDING; 188 PetscCall(VecCopy(th->X1, ts->vec_sol)); 189 PetscCall(TSAdaptChoose(ts->adapt, ts, ts->time_step, NULL, &next_time_step, &accept)); 190 th->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 191 if (!accept) { 192 PetscCall(VecCopy(th->X0, ts->vec_sol)); 193 ts->time_step = next_time_step; 194 goto reject_step; 195 } 196 197 ts->ptime += ts->time_step; 198 ts->time_step = next_time_step; 199 break; 200 201 reject_step: 202 ts->reject++; 203 accept = PETSC_FALSE; 204 if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 205 ts->reason = TS_DIVERGED_STEP_REJECTED; 206 PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, rejections)); 207 } 208 } 209 PetscFunctionReturn(0); 210 } 211 212 static PetscErrorCode TSEvaluateWLTE_Alpha(TS ts, NormType wnormtype, PetscInt *order, PetscReal *wlte) { 213 TS_Alpha *th = (TS_Alpha *)ts->data; 214 Vec X = th->X1; /* X = solution */ 215 Vec Y = th->vec_lte_work; /* Y = X + LTE */ 216 PetscReal wltea, wlter; 217 218 PetscFunctionBegin; 219 if (!th->vec_sol_prev) { 220 *wlte = -1; 221 PetscFunctionReturn(0); 222 } 223 if (!th->vec_lte_work) { 224 *wlte = -1; 225 PetscFunctionReturn(0); 226 } 227 if (ts->steprestart) { 228 /* th->vec_lte_work is set to the LTE in TSAlpha_Restart() */ 229 PetscCall(VecAXPY(Y, 1, X)); 230 } else { 231 /* Compute LTE using backward differences with non-constant time step */ 232 PetscReal h = ts->time_step, h_prev = ts->ptime - ts->ptime_prev; 233 PetscReal a = 1 + h_prev / h; 234 PetscScalar scal[3]; 235 Vec vecs[3]; 236 scal[0] = +1 / a; 237 scal[1] = -1 / (a - 1); 238 scal[2] = +1 / (a * (a - 1)); 239 vecs[0] = th->X1; 240 vecs[1] = th->X0; 241 vecs[2] = th->vec_sol_prev; 242 PetscCall(VecCopy(X, Y)); 243 PetscCall(VecMAXPY(Y, 3, scal, vecs)); 244 } 245 PetscCall(TSErrorWeightedNorm(ts, X, Y, wnormtype, wlte, &wltea, &wlter)); 246 if (order) *order = 2; 247 PetscFunctionReturn(0); 248 } 249 250 static PetscErrorCode TSRollBack_Alpha(TS ts) { 251 TS_Alpha *th = (TS_Alpha *)ts->data; 252 253 PetscFunctionBegin; 254 PetscCall(VecCopy(th->X0, ts->vec_sol)); 255 PetscFunctionReturn(0); 256 } 257 258 static PetscErrorCode TSInterpolate_Alpha(TS ts, PetscReal t, Vec X) { 259 TS_Alpha *th = (TS_Alpha *)ts->data; 260 PetscReal dt = t - ts->ptime; 261 262 PetscFunctionBegin; 263 PetscCall(VecCopy(ts->vec_sol, X)); 264 PetscCall(VecAXPY(X, th->Gamma * dt, th->V1)); 265 PetscCall(VecAXPY(X, (1 - th->Gamma) * dt, th->V0)); 266 PetscFunctionReturn(0); 267 } 268 269 static PetscErrorCode SNESTSFormFunction_Alpha(PETSC_UNUSED SNES snes, Vec X, Vec F, TS ts) { 270 TS_Alpha *th = (TS_Alpha *)ts->data; 271 PetscReal ta = th->stage_time; 272 Vec Xa = th->Xa, Va = th->Va; 273 274 PetscFunctionBegin; 275 PetscCall(TSAlpha_StageVecs(ts, X)); 276 /* F = Function(ta,Xa,Va) */ 277 PetscCall(TSComputeIFunction(ts, ta, Xa, Va, F, PETSC_FALSE)); 278 PetscCall(VecScale(F, th->scale_F)); 279 PetscFunctionReturn(0); 280 } 281 282 static PetscErrorCode SNESTSFormJacobian_Alpha(PETSC_UNUSED SNES snes, PETSC_UNUSED Vec X, Mat J, Mat P, TS ts) { 283 TS_Alpha *th = (TS_Alpha *)ts->data; 284 PetscReal ta = th->stage_time; 285 Vec Xa = th->Xa, Va = th->Va; 286 PetscReal dVdX = th->shift_V; 287 288 PetscFunctionBegin; 289 /* J,P = Jacobian(ta,Xa,Va) */ 290 PetscCall(TSComputeIJacobian(ts, ta, Xa, Va, dVdX, J, P, PETSC_FALSE)); 291 PetscFunctionReturn(0); 292 } 293 294 static PetscErrorCode TSReset_Alpha(TS ts) { 295 TS_Alpha *th = (TS_Alpha *)ts->data; 296 297 PetscFunctionBegin; 298 PetscCall(VecDestroy(&th->X0)); 299 PetscCall(VecDestroy(&th->Xa)); 300 PetscCall(VecDestroy(&th->X1)); 301 PetscCall(VecDestroy(&th->V0)); 302 PetscCall(VecDestroy(&th->Va)); 303 PetscCall(VecDestroy(&th->V1)); 304 PetscCall(VecDestroy(&th->vec_sol_prev)); 305 PetscCall(VecDestroy(&th->vec_lte_work)); 306 PetscFunctionReturn(0); 307 } 308 309 static PetscErrorCode TSDestroy_Alpha(TS ts) { 310 PetscFunctionBegin; 311 PetscCall(TSReset_Alpha(ts)); 312 PetscCall(PetscFree(ts->data)); 313 314 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlphaSetRadius_C", NULL)); 315 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlphaSetParams_C", NULL)); 316 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlphaGetParams_C", NULL)); 317 PetscFunctionReturn(0); 318 } 319 320 static PetscErrorCode TSSetUp_Alpha(TS ts) { 321 TS_Alpha *th = (TS_Alpha *)ts->data; 322 PetscBool match; 323 324 PetscFunctionBegin; 325 PetscCall(VecDuplicate(ts->vec_sol, &th->X0)); 326 PetscCall(VecDuplicate(ts->vec_sol, &th->Xa)); 327 PetscCall(VecDuplicate(ts->vec_sol, &th->X1)); 328 PetscCall(VecDuplicate(ts->vec_sol, &th->V0)); 329 PetscCall(VecDuplicate(ts->vec_sol, &th->Va)); 330 PetscCall(VecDuplicate(ts->vec_sol, &th->V1)); 331 332 PetscCall(TSGetAdapt(ts, &ts->adapt)); 333 PetscCall(TSAdaptCandidatesClear(ts->adapt)); 334 PetscCall(PetscObjectTypeCompare((PetscObject)ts->adapt, TSADAPTNONE, &match)); 335 if (!match) { 336 PetscCall(VecDuplicate(ts->vec_sol, &th->vec_sol_prev)); 337 PetscCall(VecDuplicate(ts->vec_sol, &th->vec_lte_work)); 338 } 339 340 PetscCall(TSGetSNES(ts, &ts->snes)); 341 PetscFunctionReturn(0); 342 } 343 344 static PetscErrorCode TSSetFromOptions_Alpha(TS ts, PetscOptionItems *PetscOptionsObject) { 345 TS_Alpha *th = (TS_Alpha *)ts->data; 346 347 PetscFunctionBegin; 348 PetscOptionsHeadBegin(PetscOptionsObject, "Generalized-Alpha ODE solver options"); 349 { 350 PetscBool flg; 351 PetscReal radius = 1; 352 PetscCall(PetscOptionsReal("-ts_alpha_radius", "Spectral radius (high-frequency dissipation)", "TSAlphaSetRadius", radius, &radius, &flg)); 353 if (flg) PetscCall(TSAlphaSetRadius(ts, radius)); 354 PetscCall(PetscOptionsReal("-ts_alpha_alpha_m", "Algorithmic parameter alpha_m", "TSAlphaSetParams", th->Alpha_m, &th->Alpha_m, NULL)); 355 PetscCall(PetscOptionsReal("-ts_alpha_alpha_f", "Algorithmic parameter alpha_f", "TSAlphaSetParams", th->Alpha_f, &th->Alpha_f, NULL)); 356 PetscCall(PetscOptionsReal("-ts_alpha_gamma", "Algorithmic parameter gamma", "TSAlphaSetParams", th->Gamma, &th->Gamma, NULL)); 357 PetscCall(TSAlphaSetParams(ts, th->Alpha_m, th->Alpha_f, th->Gamma)); 358 } 359 PetscOptionsHeadEnd(); 360 PetscFunctionReturn(0); 361 } 362 363 static PetscErrorCode TSView_Alpha(TS ts, PetscViewer viewer) { 364 TS_Alpha *th = (TS_Alpha *)ts->data; 365 PetscBool iascii; 366 367 PetscFunctionBegin; 368 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 369 if (iascii) PetscCall(PetscViewerASCIIPrintf(viewer, " Alpha_m=%g, Alpha_f=%g, Gamma=%g\n", (double)th->Alpha_m, (double)th->Alpha_f, (double)th->Gamma)); 370 PetscFunctionReturn(0); 371 } 372 373 static PetscErrorCode TSAlphaSetRadius_Alpha(TS ts, PetscReal radius) { 374 PetscReal alpha_m, alpha_f, gamma; 375 376 PetscFunctionBegin; 377 PetscCheck(radius >= 0 && radius <= 1, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_OUTOFRANGE, "Radius %g not in range [0,1]", (double)radius); 378 alpha_m = (PetscReal)0.5 * (3 - radius) / (1 + radius); 379 alpha_f = 1 / (1 + radius); 380 gamma = (PetscReal)0.5 + alpha_m - alpha_f; 381 PetscCall(TSAlphaSetParams(ts, alpha_m, alpha_f, gamma)); 382 PetscFunctionReturn(0); 383 } 384 385 static PetscErrorCode TSAlphaSetParams_Alpha(TS ts, PetscReal alpha_m, PetscReal alpha_f, PetscReal gamma) { 386 TS_Alpha *th = (TS_Alpha *)ts->data; 387 PetscReal tol = 100 * PETSC_MACHINE_EPSILON; 388 PetscReal res = ((PetscReal)0.5 + alpha_m - alpha_f) - gamma; 389 390 PetscFunctionBegin; 391 th->Alpha_m = alpha_m; 392 th->Alpha_f = alpha_f; 393 th->Gamma = gamma; 394 th->order = (PetscAbsReal(res) < tol) ? 2 : 1; 395 PetscFunctionReturn(0); 396 } 397 398 static PetscErrorCode TSAlphaGetParams_Alpha(TS ts, PetscReal *alpha_m, PetscReal *alpha_f, PetscReal *gamma) { 399 TS_Alpha *th = (TS_Alpha *)ts->data; 400 401 PetscFunctionBegin; 402 if (alpha_m) *alpha_m = th->Alpha_m; 403 if (alpha_f) *alpha_f = th->Alpha_f; 404 if (gamma) *gamma = th->Gamma; 405 PetscFunctionReturn(0); 406 } 407 408 /*MC 409 TSALPHA - ODE/DAE solver using the implicit Generalized-Alpha method 410 for first-order systems 411 412 Level: beginner 413 414 References: 415 + * - K.E. Jansen, C.H. Whiting, G.M. Hulber, "A generalized-alpha 416 method for integrating the filtered Navier-Stokes equations with a 417 stabilized finite element method", Computer Methods in Applied 418 Mechanics and Engineering, 190, 305-319, 2000. 419 DOI: 10.1016/S0045-7825(00)00203-6. 420 - * - J. Chung, G.M.Hubert. "A Time Integration Algorithm for Structural 421 Dynamics with Improved Numerical Dissipation: The Generalized-alpha 422 Method" ASME Journal of Applied Mechanics, 60, 371:375, 1993. 423 424 .seealso: `TS`, `TSCreate()`, `TSSetType()`, `TSAlphaSetRadius()`, `TSAlphaSetParams()` 425 M*/ 426 PETSC_EXTERN PetscErrorCode TSCreate_Alpha(TS ts) { 427 TS_Alpha *th; 428 429 PetscFunctionBegin; 430 ts->ops->reset = TSReset_Alpha; 431 ts->ops->destroy = TSDestroy_Alpha; 432 ts->ops->view = TSView_Alpha; 433 ts->ops->setup = TSSetUp_Alpha; 434 ts->ops->setfromoptions = TSSetFromOptions_Alpha; 435 ts->ops->step = TSStep_Alpha; 436 ts->ops->evaluatewlte = TSEvaluateWLTE_Alpha; 437 ts->ops->rollback = TSRollBack_Alpha; 438 ts->ops->interpolate = TSInterpolate_Alpha; 439 ts->ops->snesfunction = SNESTSFormFunction_Alpha; 440 ts->ops->snesjacobian = SNESTSFormJacobian_Alpha; 441 ts->default_adapt_type = TSADAPTNONE; 442 443 ts->usessnes = PETSC_TRUE; 444 445 PetscCall(PetscNewLog(ts, &th)); 446 ts->data = (void *)th; 447 448 th->Alpha_m = 0.5; 449 th->Alpha_f = 0.5; 450 th->Gamma = 0.5; 451 th->order = 2; 452 453 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlphaSetRadius_C", TSAlphaSetRadius_Alpha)); 454 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlphaSetParams_C", TSAlphaSetParams_Alpha)); 455 PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlphaGetParams_C", TSAlphaGetParams_Alpha)); 456 PetscFunctionReturn(0); 457 } 458 459 /*@ 460 TSAlphaSetRadius - sets the desired spectral radius of the method 461 (i.e. high-frequency numerical damping) 462 463 Logically Collective on TS 464 465 The algorithmic parameters \alpha_m and \alpha_f of the 466 generalized-\alpha method can be computed in terms of a specified 467 spectral radius \rho in [0,1] for infinite time step in order to 468 control high-frequency numerical damping: 469 \alpha_m = 0.5*(3-\rho)/(1+\rho) 470 \alpha_f = 1/(1+\rho) 471 472 Input Parameters: 473 + ts - timestepping context 474 - radius - the desired spectral radius 475 476 Options Database: 477 . -ts_alpha_radius <radius> - set alpha radius 478 479 Level: intermediate 480 481 .seealso: `TSAlphaSetParams()`, `TSAlphaGetParams()` 482 @*/ 483 PetscErrorCode TSAlphaSetRadius(TS ts, PetscReal radius) { 484 PetscFunctionBegin; 485 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 486 PetscValidLogicalCollectiveReal(ts, radius, 2); 487 PetscCheck(radius >= 0 && radius <= 1, ((PetscObject)ts)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Radius %g not in range [0,1]", (double)radius); 488 PetscTryMethod(ts, "TSAlphaSetRadius_C", (TS, PetscReal), (ts, radius)); 489 PetscFunctionReturn(0); 490 } 491 492 /*@ 493 TSAlphaSetParams - sets the algorithmic parameters for TSALPHA 494 495 Logically Collective on TS 496 497 Second-order accuracy can be obtained so long as: 498 \gamma = 0.5 + alpha_m - alpha_f 499 500 Unconditional stability requires: 501 \alpha_m >= \alpha_f >= 0.5 502 503 Backward Euler method is recovered with: 504 \alpha_m = \alpha_f = gamma = 1 505 506 Input Parameters: 507 + ts - timestepping context 508 . \alpha_m - algorithmic parameter 509 . \alpha_f - algorithmic parameter 510 - \gamma - algorithmic parameter 511 512 Options Database: 513 + -ts_alpha_alpha_m <alpha_m> - set alpha_m 514 . -ts_alpha_alpha_f <alpha_f> - set alpha_f 515 - -ts_alpha_gamma <gamma> - set gamma 516 517 Note: 518 Use of this function is normally only required to hack TSALPHA to 519 use a modified integration scheme. Users should call 520 TSAlphaSetRadius() to set the desired spectral radius of the methods 521 (i.e. high-frequency damping) in order so select optimal values for 522 these parameters. 523 524 Level: advanced 525 526 .seealso: `TSAlphaSetRadius()`, `TSAlphaGetParams()` 527 @*/ 528 PetscErrorCode TSAlphaSetParams(TS ts, PetscReal alpha_m, PetscReal alpha_f, PetscReal gamma) { 529 PetscFunctionBegin; 530 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 531 PetscValidLogicalCollectiveReal(ts, alpha_m, 2); 532 PetscValidLogicalCollectiveReal(ts, alpha_f, 3); 533 PetscValidLogicalCollectiveReal(ts, gamma, 4); 534 PetscTryMethod(ts, "TSAlphaSetParams_C", (TS, PetscReal, PetscReal, PetscReal), (ts, alpha_m, alpha_f, gamma)); 535 PetscFunctionReturn(0); 536 } 537 538 /*@ 539 TSAlphaGetParams - gets the algorithmic parameters for TSALPHA 540 541 Not Collective 542 543 Input Parameter: 544 . ts - timestepping context 545 546 Output Parameters: 547 + \alpha_m - algorithmic parameter 548 . \alpha_f - algorithmic parameter 549 - \gamma - algorithmic parameter 550 551 Note: 552 Use of this function is normally only required to hack TSALPHA to 553 use a modified integration scheme. Users should call 554 TSAlphaSetRadius() to set the high-frequency damping (i.e. spectral 555 radius of the method) in order so select optimal values for these 556 parameters. 557 558 Level: advanced 559 560 .seealso: `TSAlphaSetRadius()`, `TSAlphaSetParams()` 561 @*/ 562 PetscErrorCode TSAlphaGetParams(TS ts, PetscReal *alpha_m, PetscReal *alpha_f, PetscReal *gamma) { 563 PetscFunctionBegin; 564 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 565 if (alpha_m) PetscValidRealPointer(alpha_m, 2); 566 if (alpha_f) PetscValidRealPointer(alpha_f, 3); 567 if (gamma) PetscValidRealPointer(gamma, 4); 568 PetscUseMethod(ts, "TSAlphaGetParams_C", (TS, PetscReal *, PetscReal *, PetscReal *), (ts, alpha_m, alpha_f, gamma)); 569 PetscFunctionReturn(0); 570 } 571