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