1818efac9SLisandro Dalcin /* 2818efac9SLisandro Dalcin Code for timestepping with implicit generalized-\alpha method 3818efac9SLisandro Dalcin for second order systems. 4818efac9SLisandro Dalcin */ 5818efac9SLisandro Dalcin #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 6818efac9SLisandro Dalcin 7818efac9SLisandro Dalcin static PetscBool cited = PETSC_FALSE; 89371c9d4SSatish Balay static const char citation[] = "@article{Chung1993,\n" 9818efac9SLisandro Dalcin " title = {A Time Integration Algorithm for Structural Dynamics with Improved Numerical Dissipation: The Generalized-$\\alpha$ Method},\n" 10818efac9SLisandro Dalcin " author = {J. Chung, G. M. Hubert},\n" 11818efac9SLisandro Dalcin " journal = {ASME Journal of Applied Mechanics},\n" 12818efac9SLisandro Dalcin " volume = {60},\n" 13818efac9SLisandro Dalcin " number = {2},\n" 14818efac9SLisandro Dalcin " pages = {371--375},\n" 15818efac9SLisandro Dalcin " year = {1993},\n" 16818efac9SLisandro Dalcin " issn = {0021-8936},\n" 17818efac9SLisandro Dalcin " doi = {http://dx.doi.org/10.1115/1.2900803}\n}\n"; 18818efac9SLisandro Dalcin 19818efac9SLisandro Dalcin typedef struct { 20818efac9SLisandro Dalcin PetscReal stage_time; 21818efac9SLisandro Dalcin PetscReal shift_V; 22818efac9SLisandro Dalcin PetscReal shift_A; 23818efac9SLisandro Dalcin PetscReal scale_F; 24818efac9SLisandro Dalcin Vec X0, Xa, X1; 25818efac9SLisandro Dalcin Vec V0, Va, V1; 26818efac9SLisandro Dalcin Vec A0, Aa, A1; 27818efac9SLisandro Dalcin 28818efac9SLisandro Dalcin Vec vec_dot; 29818efac9SLisandro Dalcin 30818efac9SLisandro Dalcin PetscReal Alpha_m; 31818efac9SLisandro Dalcin PetscReal Alpha_f; 32818efac9SLisandro Dalcin PetscReal Gamma; 33818efac9SLisandro Dalcin PetscReal Beta; 34818efac9SLisandro Dalcin PetscInt order; 35818efac9SLisandro Dalcin 36818efac9SLisandro Dalcin Vec vec_sol_prev; 37818efac9SLisandro Dalcin Vec vec_dot_prev; 38818efac9SLisandro Dalcin Vec vec_lte_work[2]; 39818efac9SLisandro Dalcin 40818efac9SLisandro Dalcin TSStepStatus status; 41*220f924aSDavid Kamensky 42*220f924aSDavid Kamensky TSAlpha2Predictor predictor; 43*220f924aSDavid Kamensky void *predictor_ctx; 44818efac9SLisandro Dalcin } TS_Alpha; 45818efac9SLisandro Dalcin 46*220f924aSDavid Kamensky /*@C 47*220f924aSDavid Kamensky TSAlpha2SetPredictor - sets the callback for computing a predictor (i.e., initial guess 48*220f924aSDavid Kamensky for the nonlinear solver). 49*220f924aSDavid Kamensky 50*220f924aSDavid Kamensky Input Parameters: 51*220f924aSDavid Kamensky + ts - timestepping context 52*220f924aSDavid Kamensky . predictor - callback to set the predictor in each step 53*220f924aSDavid Kamensky - ctx - the application context, which may be set to `NULL` if not used 54*220f924aSDavid Kamensky 55*220f924aSDavid Kamensky Level: intermediate 56*220f924aSDavid Kamensky 57*220f924aSDavid Kamensky Notes: 58*220f924aSDavid Kamensky 59*220f924aSDavid Kamensky If this function is never called, a same-state-vector predictor will be used, i.e., 60*220f924aSDavid Kamensky the initial guess will be the converged solution from the previous time step, without regard 61*220f924aSDavid Kamensky for the previous velocity or acceleration. 62*220f924aSDavid Kamensky 63*220f924aSDavid Kamensky .seealso: [](ch_ts), `TS`, `TSALPHA2`, `TSAlpha2Predictor` 64*220f924aSDavid Kamensky @*/ 65*220f924aSDavid Kamensky PetscErrorCode TSAlpha2SetPredictor(TS ts, TSAlpha2Predictor predictor, void *ctx) 66*220f924aSDavid Kamensky { 67*220f924aSDavid Kamensky TS_Alpha *th = (TS_Alpha *)(ts->data); 68*220f924aSDavid Kamensky 69*220f924aSDavid Kamensky PetscFunctionBegin; 70*220f924aSDavid Kamensky th->predictor = predictor; 71*220f924aSDavid Kamensky th->predictor_ctx = ctx; 72*220f924aSDavid Kamensky PetscFunctionReturn(PETSC_SUCCESS); 73*220f924aSDavid Kamensky } 74*220f924aSDavid Kamensky 75*220f924aSDavid Kamensky static PetscErrorCode TSAlpha_ApplyPredictor(TS ts, Vec X1) 76*220f924aSDavid Kamensky { 77*220f924aSDavid Kamensky /* Apply a custom predictor if set, or default to same-displacement. */ 78*220f924aSDavid Kamensky TS_Alpha *th = (TS_Alpha *)(ts->data); 79*220f924aSDavid Kamensky 80*220f924aSDavid Kamensky PetscFunctionBegin; 81*220f924aSDavid Kamensky if (th->predictor) PetscCall(th->predictor(ts, th->X0, th->V0, th->A0, X1, th->predictor_ctx)); 82*220f924aSDavid Kamensky else PetscCall(VecCopy(th->X0, X1)); 83*220f924aSDavid Kamensky PetscFunctionReturn(PETSC_SUCCESS); 84*220f924aSDavid Kamensky } 85*220f924aSDavid Kamensky 86d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAlpha_StageTime(TS ts) 87d71ae5a4SJacob Faibussowitsch { 88818efac9SLisandro Dalcin TS_Alpha *th = (TS_Alpha *)ts->data; 89818efac9SLisandro Dalcin PetscReal t = ts->ptime; 90818efac9SLisandro Dalcin PetscReal dt = ts->time_step; 91818efac9SLisandro Dalcin PetscReal Alpha_m = th->Alpha_m; 92818efac9SLisandro Dalcin PetscReal Alpha_f = th->Alpha_f; 93818efac9SLisandro Dalcin PetscReal Gamma = th->Gamma; 94818efac9SLisandro Dalcin PetscReal Beta = th->Beta; 95818efac9SLisandro Dalcin 96818efac9SLisandro Dalcin PetscFunctionBegin; 97818efac9SLisandro Dalcin th->stage_time = t + Alpha_f * dt; 98818efac9SLisandro Dalcin th->shift_V = Gamma / (dt * Beta); 99818efac9SLisandro Dalcin th->shift_A = Alpha_m / (Alpha_f * dt * dt * Beta); 100818efac9SLisandro Dalcin th->scale_F = 1 / Alpha_f; 1013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 102818efac9SLisandro Dalcin } 103818efac9SLisandro Dalcin 104d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAlpha_StageVecs(TS ts, Vec X) 105d71ae5a4SJacob Faibussowitsch { 106818efac9SLisandro Dalcin TS_Alpha *th = (TS_Alpha *)ts->data; 107818efac9SLisandro Dalcin Vec X1 = X, V1 = th->V1, A1 = th->A1; 108818efac9SLisandro Dalcin Vec Xa = th->Xa, Va = th->Va, Aa = th->Aa; 109818efac9SLisandro Dalcin Vec X0 = th->X0, V0 = th->V0, A0 = th->A0; 110818efac9SLisandro Dalcin PetscReal dt = ts->time_step; 111818efac9SLisandro Dalcin PetscReal Alpha_m = th->Alpha_m; 112818efac9SLisandro Dalcin PetscReal Alpha_f = th->Alpha_f; 113818efac9SLisandro Dalcin PetscReal Gamma = th->Gamma; 114818efac9SLisandro Dalcin PetscReal Beta = th->Beta; 115818efac9SLisandro Dalcin 116818efac9SLisandro Dalcin PetscFunctionBegin; 117818efac9SLisandro Dalcin /* A1 = ... */ 1189566063dSJacob Faibussowitsch PetscCall(VecWAXPY(A1, -1.0, X0, X1)); 1199566063dSJacob Faibussowitsch PetscCall(VecAXPY(A1, -dt, V0)); 1209566063dSJacob Faibussowitsch PetscCall(VecAXPBY(A1, -(1 - 2 * Beta) / (2 * Beta), 1 / (dt * dt * Beta), A0)); 121818efac9SLisandro Dalcin /* V1 = ... */ 1229566063dSJacob Faibussowitsch PetscCall(VecWAXPY(V1, (1.0 - Gamma) / Gamma, A0, A1)); 1239566063dSJacob Faibussowitsch PetscCall(VecAYPX(V1, dt * Gamma, V0)); 124818efac9SLisandro Dalcin /* Xa = X0 + Alpha_f*(X1-X0) */ 1259566063dSJacob Faibussowitsch PetscCall(VecWAXPY(Xa, -1.0, X0, X1)); 1269566063dSJacob Faibussowitsch PetscCall(VecAYPX(Xa, Alpha_f, X0)); 127818efac9SLisandro Dalcin /* Va = V0 + Alpha_f*(V1-V0) */ 1289566063dSJacob Faibussowitsch PetscCall(VecWAXPY(Va, -1.0, V0, V1)); 1299566063dSJacob Faibussowitsch PetscCall(VecAYPX(Va, Alpha_f, V0)); 130818efac9SLisandro Dalcin /* Aa = A0 + Alpha_m*(A1-A0) */ 1319566063dSJacob Faibussowitsch PetscCall(VecWAXPY(Aa, -1.0, A0, A1)); 1329566063dSJacob Faibussowitsch PetscCall(VecAYPX(Aa, Alpha_m, A0)); 1333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 134818efac9SLisandro Dalcin } 135818efac9SLisandro Dalcin 136d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAlpha_SNESSolve(TS ts, Vec b, Vec x) 137d71ae5a4SJacob Faibussowitsch { 138818efac9SLisandro Dalcin PetscInt nits, lits; 139818efac9SLisandro Dalcin 140818efac9SLisandro Dalcin PetscFunctionBegin; 1419566063dSJacob Faibussowitsch PetscCall(SNESSolve(ts->snes, b, x)); 1429566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(ts->snes, &nits)); 1439566063dSJacob Faibussowitsch PetscCall(SNESGetLinearSolveIterations(ts->snes, &lits)); 1449371c9d4SSatish Balay ts->snes_its += nits; 1459371c9d4SSatish Balay ts->ksp_its += lits; 1463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 147818efac9SLisandro Dalcin } 148818efac9SLisandro Dalcin 149818efac9SLisandro Dalcin /* 150818efac9SLisandro Dalcin Compute a consistent initial state for the generalized-alpha method. 151818efac9SLisandro Dalcin - Solve two successive backward Euler steps with halved time step. 152818efac9SLisandro Dalcin - Compute the initial second time derivative using backward differences. 153818efac9SLisandro Dalcin - If using adaptivity, estimate the LTE of the initial step. 154818efac9SLisandro Dalcin */ 155d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAlpha_Restart(TS ts, PetscBool *initok) 156d71ae5a4SJacob Faibussowitsch { 157818efac9SLisandro Dalcin TS_Alpha *th = (TS_Alpha *)ts->data; 158818efac9SLisandro Dalcin PetscReal time_step; 159818efac9SLisandro Dalcin PetscReal alpha_m, alpha_f, gamma, beta; 160818efac9SLisandro Dalcin Vec X0 = ts->vec_sol, X1, X2 = th->X1; 161818efac9SLisandro Dalcin Vec V0 = ts->vec_dot, V1, V2 = th->V1; 162818efac9SLisandro Dalcin PetscBool stageok; 163818efac9SLisandro Dalcin 164818efac9SLisandro Dalcin PetscFunctionBegin; 1659566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X0, &X1)); 1669566063dSJacob Faibussowitsch PetscCall(VecDuplicate(V0, &V1)); 167818efac9SLisandro Dalcin 168818efac9SLisandro Dalcin /* Setup backward Euler with halved time step */ 1699566063dSJacob Faibussowitsch PetscCall(TSAlpha2GetParams(ts, &alpha_m, &alpha_f, &gamma, &beta)); 1709566063dSJacob Faibussowitsch PetscCall(TSAlpha2SetParams(ts, 1, 1, 1, 0.5)); 1719566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts, &time_step)); 172818efac9SLisandro Dalcin ts->time_step = time_step / 2; 1739566063dSJacob Faibussowitsch PetscCall(TSAlpha_StageTime(ts)); 174818efac9SLisandro Dalcin th->stage_time = ts->ptime; 1759566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(th->A0)); 176818efac9SLisandro Dalcin 177818efac9SLisandro Dalcin /* First BE step, (t0,X0,V0) -> (t1,X1,V1) */ 178818efac9SLisandro Dalcin th->stage_time += ts->time_step; 1799566063dSJacob Faibussowitsch PetscCall(VecCopy(X0, th->X0)); 1809566063dSJacob Faibussowitsch PetscCall(VecCopy(V0, th->V0)); 1819566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts, th->stage_time)); 182*220f924aSDavid Kamensky PetscCall(TSAlpha_ApplyPredictor(ts, X1)); 1839566063dSJacob Faibussowitsch PetscCall(TSAlpha_SNESSolve(ts, NULL, X1)); 1849566063dSJacob Faibussowitsch PetscCall(VecCopy(th->V1, V1)); 1859566063dSJacob Faibussowitsch PetscCall(TSPostStage(ts, th->stage_time, 0, &X1)); 1869566063dSJacob Faibussowitsch PetscCall(TSAdaptCheckStage(ts->adapt, ts, th->stage_time, X1, &stageok)); 187818efac9SLisandro Dalcin if (!stageok) goto finally; 188818efac9SLisandro Dalcin 189818efac9SLisandro Dalcin /* Second BE step, (t1,X1,V1) -> (t2,X2,V2) */ 190818efac9SLisandro Dalcin th->stage_time += ts->time_step; 1919566063dSJacob Faibussowitsch PetscCall(VecCopy(X1, th->X0)); 1929566063dSJacob Faibussowitsch PetscCall(VecCopy(V1, th->V0)); 1939566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts, th->stage_time)); 194*220f924aSDavid Kamensky PetscCall(TSAlpha_ApplyPredictor(ts, X2)); 1959566063dSJacob Faibussowitsch PetscCall(TSAlpha_SNESSolve(ts, NULL, X2)); 1969566063dSJacob Faibussowitsch PetscCall(VecCopy(th->V1, V2)); 1979566063dSJacob Faibussowitsch PetscCall(TSPostStage(ts, th->stage_time, 0, &X2)); 1989566063dSJacob Faibussowitsch PetscCall(TSAdaptCheckStage(ts->adapt, ts, th->stage_time, X1, &stageok)); 199818efac9SLisandro Dalcin if (!stageok) goto finally; 200818efac9SLisandro Dalcin 201818efac9SLisandro Dalcin /* Compute A0 ~ dV/dt at t0 with backward differences */ 2029566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(th->A0)); 2039566063dSJacob Faibussowitsch PetscCall(VecAXPY(th->A0, -3 / ts->time_step, V0)); 2049566063dSJacob Faibussowitsch PetscCall(VecAXPY(th->A0, +4 / ts->time_step, V1)); 2059566063dSJacob Faibussowitsch PetscCall(VecAXPY(th->A0, -1 / ts->time_step, V2)); 206818efac9SLisandro Dalcin 207818efac9SLisandro Dalcin /* Rough, lower-order estimate LTE of the initial step */ 2082ffb9264SLisandro Dalcin if (th->vec_lte_work[0]) { 2099566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(th->vec_lte_work[0])); 2109566063dSJacob Faibussowitsch PetscCall(VecAXPY(th->vec_lte_work[0], +2, X2)); 2119566063dSJacob Faibussowitsch PetscCall(VecAXPY(th->vec_lte_work[0], -4, X1)); 2129566063dSJacob Faibussowitsch PetscCall(VecAXPY(th->vec_lte_work[0], +2, X0)); 213818efac9SLisandro Dalcin } 2142ffb9264SLisandro Dalcin if (th->vec_lte_work[1]) { 2159566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(th->vec_lte_work[1])); 2169566063dSJacob Faibussowitsch PetscCall(VecAXPY(th->vec_lte_work[1], +2, V2)); 2179566063dSJacob Faibussowitsch PetscCall(VecAXPY(th->vec_lte_work[1], -4, V1)); 2189566063dSJacob Faibussowitsch PetscCall(VecAXPY(th->vec_lte_work[1], +2, V0)); 219818efac9SLisandro Dalcin } 220818efac9SLisandro Dalcin 221818efac9SLisandro Dalcin finally: 222818efac9SLisandro Dalcin /* Revert TSAlpha to the initial state (t0,X0,V0) */ 223818efac9SLisandro Dalcin if (initok) *initok = stageok; 2249566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, time_step)); 2259566063dSJacob Faibussowitsch PetscCall(TSAlpha2SetParams(ts, alpha_m, alpha_f, gamma, beta)); 2269566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, th->X0)); 2279566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_dot, th->V0)); 228818efac9SLisandro Dalcin 2299566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X1)); 2309566063dSJacob Faibussowitsch PetscCall(VecDestroy(&V1)); 2313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 232818efac9SLisandro Dalcin } 233818efac9SLisandro Dalcin 234d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSStep_Alpha(TS ts) 235d71ae5a4SJacob Faibussowitsch { 236818efac9SLisandro Dalcin TS_Alpha *th = (TS_Alpha *)ts->data; 237818efac9SLisandro Dalcin PetscInt rejections = 0; 238818efac9SLisandro Dalcin PetscBool stageok, accept = PETSC_TRUE; 239818efac9SLisandro Dalcin PetscReal next_time_step = ts->time_step; 240818efac9SLisandro Dalcin 241818efac9SLisandro Dalcin PetscFunctionBegin; 2429566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister(citation, &cited)); 243818efac9SLisandro Dalcin 244818efac9SLisandro Dalcin if (!ts->steprollback) { 2459566063dSJacob Faibussowitsch if (th->vec_sol_prev) PetscCall(VecCopy(th->X0, th->vec_sol_prev)); 2469566063dSJacob Faibussowitsch if (th->vec_dot_prev) PetscCall(VecCopy(th->V0, th->vec_dot_prev)); 2479566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, th->X0)); 2489566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_dot, th->V0)); 2499566063dSJacob Faibussowitsch PetscCall(VecCopy(th->A1, th->A0)); 250818efac9SLisandro Dalcin } 251818efac9SLisandro Dalcin 252818efac9SLisandro Dalcin th->status = TS_STEP_INCOMPLETE; 253818efac9SLisandro Dalcin while (!ts->reason && th->status != TS_STEP_COMPLETE) { 254818efac9SLisandro Dalcin if (ts->steprestart) { 2559566063dSJacob Faibussowitsch PetscCall(TSAlpha_Restart(ts, &stageok)); 256818efac9SLisandro Dalcin if (!stageok) goto reject_step; 257818efac9SLisandro Dalcin } 258818efac9SLisandro Dalcin 2599566063dSJacob Faibussowitsch PetscCall(TSAlpha_StageTime(ts)); 260*220f924aSDavid Kamensky PetscCall(TSAlpha_ApplyPredictor(ts, th->X1)); 2619566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts, th->stage_time)); 2629566063dSJacob Faibussowitsch PetscCall(TSAlpha_SNESSolve(ts, NULL, th->X1)); 2639566063dSJacob Faibussowitsch PetscCall(TSPostStage(ts, th->stage_time, 0, &th->Xa)); 2649566063dSJacob Faibussowitsch PetscCall(TSAdaptCheckStage(ts->adapt, ts, th->stage_time, th->Xa, &stageok)); 265818efac9SLisandro Dalcin if (!stageok) goto reject_step; 266818efac9SLisandro Dalcin 267818efac9SLisandro Dalcin th->status = TS_STEP_PENDING; 2689566063dSJacob Faibussowitsch PetscCall(VecCopy(th->X1, ts->vec_sol)); 2699566063dSJacob Faibussowitsch PetscCall(VecCopy(th->V1, ts->vec_dot)); 2709566063dSJacob Faibussowitsch PetscCall(TSAdaptChoose(ts->adapt, ts, ts->time_step, NULL, &next_time_step, &accept)); 271818efac9SLisandro Dalcin th->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 272818efac9SLisandro Dalcin if (!accept) { 2739566063dSJacob Faibussowitsch PetscCall(VecCopy(th->X0, ts->vec_sol)); 2749566063dSJacob Faibussowitsch PetscCall(VecCopy(th->V0, ts->vec_dot)); 275818efac9SLisandro Dalcin ts->time_step = next_time_step; 276818efac9SLisandro Dalcin goto reject_step; 277818efac9SLisandro Dalcin } 278818efac9SLisandro Dalcin 279818efac9SLisandro Dalcin ts->ptime += ts->time_step; 280818efac9SLisandro Dalcin ts->time_step = next_time_step; 281818efac9SLisandro Dalcin break; 282818efac9SLisandro Dalcin 283818efac9SLisandro Dalcin reject_step: 2849371c9d4SSatish Balay ts->reject++; 2859371c9d4SSatish Balay accept = PETSC_FALSE; 286818efac9SLisandro Dalcin if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 287818efac9SLisandro Dalcin ts->reason = TS_DIVERGED_STEP_REJECTED; 28863a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, rejections)); 289818efac9SLisandro Dalcin } 290818efac9SLisandro Dalcin } 2913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 292818efac9SLisandro Dalcin } 293818efac9SLisandro Dalcin 294d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSEvaluateWLTE_Alpha(TS ts, NormType wnormtype, PetscInt *order, PetscReal *wlte) 295d71ae5a4SJacob Faibussowitsch { 296818efac9SLisandro Dalcin TS_Alpha *th = (TS_Alpha *)ts->data; 297818efac9SLisandro Dalcin Vec X = th->X1; /* X = solution */ 298818efac9SLisandro Dalcin Vec V = th->V1; /* V = solution */ 299818efac9SLisandro Dalcin Vec Y = th->vec_lte_work[0]; /* Y = X + LTE */ 300818efac9SLisandro Dalcin Vec Z = th->vec_lte_work[1]; /* Z = V + LTE */ 3017453f775SEmil Constantinescu PetscReal enormX, enormV, enormXa, enormVa, enormXr, enormVr; 302818efac9SLisandro Dalcin 303818efac9SLisandro Dalcin PetscFunctionBegin; 3049371c9d4SSatish Balay if (!th->vec_sol_prev) { 3059371c9d4SSatish Balay *wlte = -1; 3063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3079371c9d4SSatish Balay } 3089371c9d4SSatish Balay if (!th->vec_dot_prev) { 3099371c9d4SSatish Balay *wlte = -1; 3103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3119371c9d4SSatish Balay } 3129371c9d4SSatish Balay if (!th->vec_lte_work[0]) { 3139371c9d4SSatish Balay *wlte = -1; 3143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3159371c9d4SSatish Balay } 3169371c9d4SSatish Balay if (!th->vec_lte_work[1]) { 3179371c9d4SSatish Balay *wlte = -1; 3183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3199371c9d4SSatish Balay } 320818efac9SLisandro Dalcin if (ts->steprestart) { 3212ffb9264SLisandro Dalcin /* th->vec_lte_prev is set to the LTE in TSAlpha_Restart() */ 3229566063dSJacob Faibussowitsch PetscCall(VecAXPY(Y, 1, X)); 3239566063dSJacob Faibussowitsch PetscCall(VecAXPY(Z, 1, V)); 324818efac9SLisandro Dalcin } else { 325818efac9SLisandro Dalcin /* Compute LTE using backward differences with non-constant time step */ 326818efac9SLisandro Dalcin PetscReal h = ts->time_step, h_prev = ts->ptime - ts->ptime_prev; 327818efac9SLisandro Dalcin PetscReal a = 1 + h_prev / h; 3289371c9d4SSatish Balay PetscScalar scal[3]; 3299371c9d4SSatish Balay Vec vecX[3], vecV[3]; 3309371c9d4SSatish Balay scal[0] = +1 / a; 3319371c9d4SSatish Balay scal[1] = -1 / (a - 1); 3329371c9d4SSatish Balay scal[2] = +1 / (a * (a - 1)); 3339371c9d4SSatish Balay vecX[0] = th->X1; 3349371c9d4SSatish Balay vecX[1] = th->X0; 3359371c9d4SSatish Balay vecX[2] = th->vec_sol_prev; 3369371c9d4SSatish Balay vecV[0] = th->V1; 3379371c9d4SSatish Balay vecV[1] = th->V0; 3389371c9d4SSatish Balay vecV[2] = th->vec_dot_prev; 3399566063dSJacob Faibussowitsch PetscCall(VecCopy(X, Y)); 3409566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Y, 3, scal, vecX)); 3419566063dSJacob Faibussowitsch PetscCall(VecCopy(V, Z)); 3429566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Z, 3, scal, vecV)); 343818efac9SLisandro Dalcin } 344818efac9SLisandro Dalcin /* XXX ts->atol and ts->vatol are not appropriate for computing enormV */ 3459566063dSJacob Faibussowitsch PetscCall(TSErrorWeightedNorm(ts, X, Y, wnormtype, &enormX, &enormXa, &enormXr)); 3469566063dSJacob Faibussowitsch PetscCall(TSErrorWeightedNorm(ts, V, Z, wnormtype, &enormV, &enormVa, &enormVr)); 3479371c9d4SSatish Balay if (wnormtype == NORM_2) *wlte = PetscSqrtReal(PetscSqr(enormX) / 2 + PetscSqr(enormV) / 2); 3489371c9d4SSatish Balay else *wlte = PetscMax(enormX, enormV); 349818efac9SLisandro Dalcin if (order) *order = 2; 3503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 351818efac9SLisandro Dalcin } 352818efac9SLisandro Dalcin 353d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRollBack_Alpha(TS ts) 354d71ae5a4SJacob Faibussowitsch { 355818efac9SLisandro Dalcin TS_Alpha *th = (TS_Alpha *)ts->data; 356818efac9SLisandro Dalcin 357818efac9SLisandro Dalcin PetscFunctionBegin; 3589566063dSJacob Faibussowitsch PetscCall(VecCopy(th->X0, ts->vec_sol)); 3599566063dSJacob Faibussowitsch PetscCall(VecCopy(th->V0, ts->vec_dot)); 3603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 361818efac9SLisandro Dalcin } 362818efac9SLisandro Dalcin 363818efac9SLisandro Dalcin /* 364818efac9SLisandro Dalcin static PetscErrorCode TSInterpolate_Alpha(TS ts,PetscReal t,Vec X,Vec V) 365818efac9SLisandro Dalcin { 366818efac9SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 367818efac9SLisandro Dalcin PetscReal dt = t - ts->ptime; 368818efac9SLisandro Dalcin 369818efac9SLisandro Dalcin PetscFunctionBegin; 3709566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_dot,V)); 3719566063dSJacob Faibussowitsch PetscCall(VecAXPY(V,dt*(1-th->Gamma),th->A0)); 3729566063dSJacob Faibussowitsch PetscCall(VecAXPY(V,dt*th->Gamma,th->A1)); 3739566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol,X)); 3749566063dSJacob Faibussowitsch PetscCall(VecAXPY(X,dt,V)); 3759566063dSJacob Faibussowitsch PetscCall(VecAXPY(X,dt*dt*((PetscReal)0.5-th->Beta),th->A0)); 3769566063dSJacob Faibussowitsch PetscCall(VecAXPY(X,dt*dt*th->Beta,th->A1)); 3773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 378818efac9SLisandro Dalcin } 379818efac9SLisandro Dalcin */ 380818efac9SLisandro Dalcin 381d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormFunction_Alpha(PETSC_UNUSED SNES snes, Vec X, Vec F, TS ts) 382d71ae5a4SJacob Faibussowitsch { 383818efac9SLisandro Dalcin TS_Alpha *th = (TS_Alpha *)ts->data; 384818efac9SLisandro Dalcin PetscReal ta = th->stage_time; 385818efac9SLisandro Dalcin Vec Xa = th->Xa, Va = th->Va, Aa = th->Aa; 386818efac9SLisandro Dalcin 387818efac9SLisandro Dalcin PetscFunctionBegin; 3889566063dSJacob Faibussowitsch PetscCall(TSAlpha_StageVecs(ts, X)); 389818efac9SLisandro Dalcin /* F = Function(ta,Xa,Va,Aa) */ 3909566063dSJacob Faibussowitsch PetscCall(TSComputeI2Function(ts, ta, Xa, Va, Aa, F)); 3919566063dSJacob Faibussowitsch PetscCall(VecScale(F, th->scale_F)); 3923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 393818efac9SLisandro Dalcin } 394818efac9SLisandro Dalcin 395d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormJacobian_Alpha(PETSC_UNUSED SNES snes, PETSC_UNUSED Vec X, Mat J, Mat P, TS ts) 396d71ae5a4SJacob Faibussowitsch { 397818efac9SLisandro Dalcin TS_Alpha *th = (TS_Alpha *)ts->data; 398818efac9SLisandro Dalcin PetscReal ta = th->stage_time; 399818efac9SLisandro Dalcin Vec Xa = th->Xa, Va = th->Va, Aa = th->Aa; 400818efac9SLisandro Dalcin PetscReal dVdX = th->shift_V, dAdX = th->shift_A; 401818efac9SLisandro Dalcin 402818efac9SLisandro Dalcin PetscFunctionBegin; 403818efac9SLisandro Dalcin /* J,P = Jacobian(ta,Xa,Va,Aa) */ 4049566063dSJacob Faibussowitsch PetscCall(TSComputeI2Jacobian(ts, ta, Xa, Va, Aa, dVdX, dAdX, J, P)); 4053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 406818efac9SLisandro Dalcin } 407818efac9SLisandro Dalcin 408d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSReset_Alpha(TS ts) 409d71ae5a4SJacob Faibussowitsch { 410818efac9SLisandro Dalcin TS_Alpha *th = (TS_Alpha *)ts->data; 411818efac9SLisandro Dalcin 412818efac9SLisandro Dalcin PetscFunctionBegin; 4139566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->X0)); 4149566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->Xa)); 4159566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->X1)); 4169566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->V0)); 4179566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->Va)); 4189566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->V1)); 4199566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->A0)); 4209566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->Aa)); 4219566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->A1)); 4229566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->vec_sol_prev)); 4239566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->vec_dot_prev)); 4249566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->vec_lte_work[0])); 4259566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->vec_lte_work[1])); 4263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 427818efac9SLisandro Dalcin } 428818efac9SLisandro Dalcin 429d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSDestroy_Alpha(TS ts) 430d71ae5a4SJacob Faibussowitsch { 431818efac9SLisandro Dalcin PetscFunctionBegin; 4329566063dSJacob Faibussowitsch PetscCall(TSReset_Alpha(ts)); 4339566063dSJacob Faibussowitsch PetscCall(PetscFree(ts->data)); 434818efac9SLisandro Dalcin 4359566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlpha2SetRadius_C", NULL)); 4369566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlpha2SetParams_C", NULL)); 4379566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlpha2GetParams_C", NULL)); 4383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 439818efac9SLisandro Dalcin } 440818efac9SLisandro Dalcin 441d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_Alpha(TS ts) 442d71ae5a4SJacob Faibussowitsch { 443818efac9SLisandro Dalcin TS_Alpha *th = (TS_Alpha *)ts->data; 4442ffb9264SLisandro Dalcin PetscBool match; 445818efac9SLisandro Dalcin 446818efac9SLisandro Dalcin PetscFunctionBegin; 4479566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &th->X0)); 4489566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &th->Xa)); 4499566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &th->X1)); 4509566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &th->V0)); 4519566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &th->Va)); 4529566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &th->V1)); 4539566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &th->A0)); 4549566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &th->Aa)); 4559566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &th->A1)); 456818efac9SLisandro Dalcin 4579566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts, &ts->adapt)); 4589566063dSJacob Faibussowitsch PetscCall(TSAdaptCandidatesClear(ts->adapt)); 4599566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)ts->adapt, TSADAPTNONE, &match)); 4602ffb9264SLisandro Dalcin if (!match) { 4619566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &th->vec_sol_prev)); 4629566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &th->vec_dot_prev)); 4639566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &th->vec_lte_work[0])); 4649566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &th->vec_lte_work[1])); 465818efac9SLisandro Dalcin } 466818efac9SLisandro Dalcin 4679566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &ts->snes)); 4683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 469818efac9SLisandro Dalcin } 470818efac9SLisandro Dalcin 471d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetFromOptions_Alpha(TS ts, PetscOptionItems *PetscOptionsObject) 472d71ae5a4SJacob Faibussowitsch { 473818efac9SLisandro Dalcin TS_Alpha *th = (TS_Alpha *)ts->data; 474818efac9SLisandro Dalcin 475818efac9SLisandro Dalcin PetscFunctionBegin; 476d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Generalized-Alpha ODE solver options"); 477818efac9SLisandro Dalcin { 478818efac9SLisandro Dalcin PetscBool flg; 479818efac9SLisandro Dalcin PetscReal radius = 1; 4809566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_alpha_radius", "Spectral radius (high-frequency dissipation)", "TSAlpha2SetRadius", radius, &radius, &flg)); 4819566063dSJacob Faibussowitsch if (flg) PetscCall(TSAlpha2SetRadius(ts, radius)); 4829566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_alpha_alpha_m", "Algorithmic parameter alpha_m", "TSAlpha2SetParams", th->Alpha_m, &th->Alpha_m, NULL)); 4839566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_alpha_alpha_f", "Algorithmic parameter alpha_f", "TSAlpha2SetParams", th->Alpha_f, &th->Alpha_f, NULL)); 4849566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_alpha_gamma", "Algorithmic parameter gamma", "TSAlpha2SetParams", th->Gamma, &th->Gamma, NULL)); 4859566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_alpha_beta", "Algorithmic parameter beta", "TSAlpha2SetParams", th->Beta, &th->Beta, NULL)); 4869566063dSJacob Faibussowitsch PetscCall(TSAlpha2SetParams(ts, th->Alpha_m, th->Alpha_f, th->Gamma, th->Beta)); 487818efac9SLisandro Dalcin } 488d0609cedSBarry Smith PetscOptionsHeadEnd(); 4893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 490818efac9SLisandro Dalcin } 491818efac9SLisandro Dalcin 492d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSView_Alpha(TS ts, PetscViewer viewer) 493d71ae5a4SJacob Faibussowitsch { 494818efac9SLisandro Dalcin TS_Alpha *th = (TS_Alpha *)ts->data; 495818efac9SLisandro Dalcin PetscBool iascii; 496818efac9SLisandro Dalcin 497818efac9SLisandro Dalcin PetscFunctionBegin; 4989566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 49948a46eb9SPierre Jolivet 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)); 5003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 501818efac9SLisandro Dalcin } 502818efac9SLisandro Dalcin 503d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAlpha2SetRadius_Alpha(TS ts, PetscReal radius) 504d71ae5a4SJacob Faibussowitsch { 505818efac9SLisandro Dalcin PetscReal alpha_m, alpha_f, gamma, beta; 506818efac9SLisandro Dalcin 507818efac9SLisandro Dalcin PetscFunctionBegin; 508cad9d221SBarry Smith PetscCheck(radius >= 0 && radius <= 1, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_OUTOFRANGE, "Radius %g not in range [0,1]", (double)radius); 509818efac9SLisandro Dalcin alpha_m = (2 - radius) / (1 + radius); 510818efac9SLisandro Dalcin alpha_f = 1 / (1 + radius); 511818efac9SLisandro Dalcin gamma = (PetscReal)0.5 + alpha_m - alpha_f; 5129371c9d4SSatish Balay beta = (PetscReal)0.5 * (1 + alpha_m - alpha_f); 5139371c9d4SSatish Balay beta *= beta; 5149566063dSJacob Faibussowitsch PetscCall(TSAlpha2SetParams(ts, alpha_m, alpha_f, gamma, beta)); 5153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 516818efac9SLisandro Dalcin } 517818efac9SLisandro Dalcin 518d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAlpha2SetParams_Alpha(TS ts, PetscReal alpha_m, PetscReal alpha_f, PetscReal gamma, PetscReal beta) 519d71ae5a4SJacob Faibussowitsch { 520818efac9SLisandro Dalcin TS_Alpha *th = (TS_Alpha *)ts->data; 521818efac9SLisandro Dalcin PetscReal tol = 100 * PETSC_MACHINE_EPSILON; 522818efac9SLisandro Dalcin PetscReal res = ((PetscReal)0.5 + alpha_m - alpha_f) - gamma; 523818efac9SLisandro Dalcin 524818efac9SLisandro Dalcin PetscFunctionBegin; 525818efac9SLisandro Dalcin th->Alpha_m = alpha_m; 526818efac9SLisandro Dalcin th->Alpha_f = alpha_f; 527818efac9SLisandro Dalcin th->Gamma = gamma; 528818efac9SLisandro Dalcin th->Beta = beta; 529818efac9SLisandro Dalcin th->order = (PetscAbsReal(res) < tol) ? 2 : 1; 5303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 531818efac9SLisandro Dalcin } 532818efac9SLisandro Dalcin 533d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAlpha2GetParams_Alpha(TS ts, PetscReal *alpha_m, PetscReal *alpha_f, PetscReal *gamma, PetscReal *beta) 534d71ae5a4SJacob Faibussowitsch { 535818efac9SLisandro Dalcin TS_Alpha *th = (TS_Alpha *)ts->data; 536818efac9SLisandro Dalcin 537818efac9SLisandro Dalcin PetscFunctionBegin; 538818efac9SLisandro Dalcin if (alpha_m) *alpha_m = th->Alpha_m; 539818efac9SLisandro Dalcin if (alpha_f) *alpha_f = th->Alpha_f; 540818efac9SLisandro Dalcin if (gamma) *gamma = th->Gamma; 541818efac9SLisandro Dalcin if (beta) *beta = th->Beta; 5423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 543818efac9SLisandro Dalcin } 544818efac9SLisandro Dalcin 545818efac9SLisandro Dalcin /*MC 54614d0ab18SJacob Faibussowitsch TSALPHA2 - ODE/DAE solver using the implicit Generalized-Alpha method for second-order systems 547818efac9SLisandro Dalcin 548818efac9SLisandro Dalcin Level: beginner 549818efac9SLisandro Dalcin 550818efac9SLisandro Dalcin References: 551606c0280SSatish Balay . * - J. Chung, G.M.Hubert. "A Time Integration Algorithm for Structural 552818efac9SLisandro Dalcin Dynamics with Improved Numerical Dissipation: The Generalized-alpha 553818efac9SLisandro Dalcin Method" ASME Journal of Applied Mechanics, 60, 371:375, 1993. 554818efac9SLisandro Dalcin 5551cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSCreate()`, `TSSetType()`, `TSAlpha2SetRadius()`, `TSAlpha2SetParams()` 556818efac9SLisandro Dalcin M*/ 557d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_Alpha2(TS ts) 558d71ae5a4SJacob Faibussowitsch { 559818efac9SLisandro Dalcin TS_Alpha *th; 560818efac9SLisandro Dalcin 561818efac9SLisandro Dalcin PetscFunctionBegin; 562818efac9SLisandro Dalcin ts->ops->reset = TSReset_Alpha; 563818efac9SLisandro Dalcin ts->ops->destroy = TSDestroy_Alpha; 564818efac9SLisandro Dalcin ts->ops->view = TSView_Alpha; 565818efac9SLisandro Dalcin ts->ops->setup = TSSetUp_Alpha; 566818efac9SLisandro Dalcin ts->ops->setfromoptions = TSSetFromOptions_Alpha; 567818efac9SLisandro Dalcin ts->ops->step = TSStep_Alpha; 568818efac9SLisandro Dalcin ts->ops->evaluatewlte = TSEvaluateWLTE_Alpha; 569818efac9SLisandro Dalcin ts->ops->rollback = TSRollBack_Alpha; 570818efac9SLisandro Dalcin /*ts->ops->interpolate = TSInterpolate_Alpha;*/ 571818efac9SLisandro Dalcin ts->ops->snesfunction = SNESTSFormFunction_Alpha; 572818efac9SLisandro Dalcin ts->ops->snesjacobian = SNESTSFormJacobian_Alpha; 5732ffb9264SLisandro Dalcin ts->default_adapt_type = TSADAPTNONE; 574818efac9SLisandro Dalcin 575efd4aadfSBarry Smith ts->usessnes = PETSC_TRUE; 576efd4aadfSBarry Smith 5774dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&th)); 578818efac9SLisandro Dalcin ts->data = (void *)th; 579818efac9SLisandro Dalcin 580818efac9SLisandro Dalcin th->Alpha_m = 0.5; 581818efac9SLisandro Dalcin th->Alpha_f = 0.5; 582818efac9SLisandro Dalcin th->Gamma = 0.5; 583818efac9SLisandro Dalcin th->Beta = 0.25; 584818efac9SLisandro Dalcin th->order = 2; 585818efac9SLisandro Dalcin 586*220f924aSDavid Kamensky th->predictor = NULL; 587*220f924aSDavid Kamensky th->predictor_ctx = NULL; 588*220f924aSDavid Kamensky 5899566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlpha2SetRadius_C", TSAlpha2SetRadius_Alpha)); 5909566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlpha2SetParams_C", TSAlpha2SetParams_Alpha)); 5919566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlpha2GetParams_C", TSAlpha2GetParams_Alpha)); 5923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 593818efac9SLisandro Dalcin } 594818efac9SLisandro Dalcin 595818efac9SLisandro Dalcin /*@ 596bcf0153eSBarry Smith TSAlpha2SetRadius - sets the desired spectral radius of the method for `TSALPHA2` 597818efac9SLisandro Dalcin (i.e. high-frequency numerical damping) 598818efac9SLisandro Dalcin 599c3339decSBarry Smith Logically Collective 600818efac9SLisandro Dalcin 601d8d19677SJose E. Roman Input Parameters: 602818efac9SLisandro Dalcin + ts - timestepping context 603818efac9SLisandro Dalcin - radius - the desired spectral radius 604818efac9SLisandro Dalcin 605bcf0153eSBarry Smith Options Database Key: 60667b8a455SSatish Balay . -ts_alpha_radius <radius> - set the desired spectral radius 607818efac9SLisandro Dalcin 608818efac9SLisandro Dalcin Level: intermediate 609818efac9SLisandro Dalcin 61014d0ab18SJacob Faibussowitsch Notes: 61114d0ab18SJacob Faibussowitsch 61214d0ab18SJacob Faibussowitsch The algorithmic parameters $\alpha_m$ and $\alpha_f$ of the generalized-$\alpha$ method can 61314d0ab18SJacob Faibussowitsch be computed in terms of a specified spectral radius $\rho$ in `[0, 1]` for infinite time step 61414d0ab18SJacob Faibussowitsch in order to control high-frequency numerical damping\: 61514d0ab18SJacob Faibussowitsch $$ 61614d0ab18SJacob Faibussowitsch \alpha_m = (2-\rho)/(1+\rho) 61714d0ab18SJacob Faibussowitsch \alpha_f = 1/(1+\rho) 61814d0ab18SJacob Faibussowitsch $$ 61914d0ab18SJacob Faibussowitsch 6201cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSALPHA2`, `TSAlpha2SetParams()`, `TSAlpha2GetParams()` 621818efac9SLisandro Dalcin @*/ 622d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAlpha2SetRadius(TS ts, PetscReal radius) 623d71ae5a4SJacob Faibussowitsch { 624818efac9SLisandro Dalcin PetscFunctionBegin; 625818efac9SLisandro Dalcin PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 626818efac9SLisandro Dalcin PetscValidLogicalCollectiveReal(ts, radius, 2); 627cad9d221SBarry Smith PetscCheck(radius >= 0 && radius <= 1, ((PetscObject)ts)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Radius %g not in range [0,1]", (double)radius); 628cac4c232SBarry Smith PetscTryMethod(ts, "TSAlpha2SetRadius_C", (TS, PetscReal), (ts, radius)); 6293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 630818efac9SLisandro Dalcin } 631818efac9SLisandro Dalcin 632818efac9SLisandro Dalcin /*@ 633bcf0153eSBarry Smith TSAlpha2SetParams - sets the algorithmic parameters for `TSALPHA2` 634818efac9SLisandro Dalcin 635c3339decSBarry Smith Logically Collective 636818efac9SLisandro Dalcin 637d8d19677SJose E. Roman Input Parameters: 638818efac9SLisandro Dalcin + ts - timestepping context 6392fe279fdSBarry Smith . alpha_m - algorithmic parameter 6402fe279fdSBarry Smith . alpha_f - algorithmic parameter 6412fe279fdSBarry Smith . gamma - algorithmic parameter 6422fe279fdSBarry Smith - beta - algorithmic parameter 643818efac9SLisandro Dalcin 644bcf0153eSBarry Smith Options Database Keys: 64567b8a455SSatish Balay + -ts_alpha_alpha_m <alpha_m> - set alpha_m 64667b8a455SSatish Balay . -ts_alpha_alpha_f <alpha_f> - set alpha_f 64767b8a455SSatish Balay . -ts_alpha_gamma <gamma> - set gamma 64867b8a455SSatish Balay - -ts_alpha_beta <beta> - set beta 649818efac9SLisandro Dalcin 650bcf0153eSBarry Smith Level: advanced 651bcf0153eSBarry Smith 65214d0ab18SJacob Faibussowitsch Notes: 65314d0ab18SJacob Faibussowitsch Second-order accuracy can be obtained so long as\: 65414d0ab18SJacob Faibussowitsch $$ 65514d0ab18SJacob Faibussowitsch \gamma = 1/2 + alpha_m - alpha_f 65614d0ab18SJacob Faibussowitsch \beta = 1/4 (1 + alpha_m - alpha_f)^2 65714d0ab18SJacob Faibussowitsch $$ 65814d0ab18SJacob Faibussowitsch 65914d0ab18SJacob Faibussowitsch Unconditional stability requires\: 66014d0ab18SJacob Faibussowitsch $$ 66114d0ab18SJacob Faibussowitsch \alpha_m >= \alpha_f >= 1/2 66214d0ab18SJacob Faibussowitsch $$ 66314d0ab18SJacob Faibussowitsch 66414d0ab18SJacob Faibussowitsch Use of this function is normally only required to hack `TSALPHA2` to use a modified 66514d0ab18SJacob Faibussowitsch integration scheme. Users should call `TSAlpha2SetRadius()` to set the desired spectral 66614d0ab18SJacob Faibussowitsch radius of the methods (i.e. high-frequency damping) in order so select optimal values for 667818efac9SLisandro Dalcin these parameters. 668818efac9SLisandro Dalcin 6691cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSALPHA2`, `TSAlpha2SetRadius()`, `TSAlpha2GetParams()` 670818efac9SLisandro Dalcin @*/ 671d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAlpha2SetParams(TS ts, PetscReal alpha_m, PetscReal alpha_f, PetscReal gamma, PetscReal beta) 672d71ae5a4SJacob Faibussowitsch { 673818efac9SLisandro Dalcin PetscFunctionBegin; 674818efac9SLisandro Dalcin PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 675818efac9SLisandro Dalcin PetscValidLogicalCollectiveReal(ts, alpha_m, 2); 676818efac9SLisandro Dalcin PetscValidLogicalCollectiveReal(ts, alpha_f, 3); 677818efac9SLisandro Dalcin PetscValidLogicalCollectiveReal(ts, gamma, 4); 678818efac9SLisandro Dalcin PetscValidLogicalCollectiveReal(ts, beta, 5); 679cac4c232SBarry Smith PetscTryMethod(ts, "TSAlpha2SetParams_C", (TS, PetscReal, PetscReal, PetscReal, PetscReal), (ts, alpha_m, alpha_f, gamma, beta)); 6803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 681818efac9SLisandro Dalcin } 682818efac9SLisandro Dalcin 683818efac9SLisandro Dalcin /*@ 684bcf0153eSBarry Smith TSAlpha2GetParams - gets the algorithmic parameters for `TSALPHA2` 685818efac9SLisandro Dalcin 686818efac9SLisandro Dalcin Not Collective 687818efac9SLisandro Dalcin 688818efac9SLisandro Dalcin Input Parameter: 689818efac9SLisandro Dalcin . ts - timestepping context 690818efac9SLisandro Dalcin 691818efac9SLisandro Dalcin Output Parameters: 6922fe279fdSBarry Smith + alpha_m - algorithmic parameter 6932fe279fdSBarry Smith . alpha_f - algorithmic parameter 6942fe279fdSBarry Smith . gamma - algorithmic parameter 6952fe279fdSBarry Smith - beta - algorithmic parameter 696818efac9SLisandro Dalcin 697bcf0153eSBarry Smith Level: advanced 698bcf0153eSBarry Smith 699818efac9SLisandro Dalcin Note: 70014d0ab18SJacob Faibussowitsch Use of this function is normally only required to hack `TSALPHA2` to use a modified 70114d0ab18SJacob Faibussowitsch integration scheme. Users should call `TSAlpha2SetRadius()` to set the high-frequency damping 70214d0ab18SJacob Faibussowitsch (i.e. spectral radius of the method) in order so select optimal values for these parameters. 703818efac9SLisandro Dalcin 7041cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSALPHA2`, `TSAlpha2SetRadius()`, `TSAlpha2SetParams()` 705818efac9SLisandro Dalcin @*/ 706d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAlpha2GetParams(TS ts, PetscReal *alpha_m, PetscReal *alpha_f, PetscReal *gamma, PetscReal *beta) 707d71ae5a4SJacob Faibussowitsch { 708818efac9SLisandro Dalcin PetscFunctionBegin; 709818efac9SLisandro Dalcin PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 7104f572ea9SToby Isaac if (alpha_m) PetscAssertPointer(alpha_m, 2); 7114f572ea9SToby Isaac if (alpha_f) PetscAssertPointer(alpha_f, 3); 7124f572ea9SToby Isaac if (gamma) PetscAssertPointer(gamma, 4); 7134f572ea9SToby Isaac if (beta) PetscAssertPointer(beta, 5); 714cac4c232SBarry Smith PetscUseMethod(ts, "TSAlpha2GetParams_C", (TS, PetscReal *, PetscReal *, PetscReal *, PetscReal *), (ts, alpha_m, alpha_f, gamma, beta)); 7153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 716818efac9SLisandro Dalcin } 717