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; 41220f924aSDavid Kamensky 428434afd1SBarry Smith TSAlpha2PredictorFn *predictor; 43220f924aSDavid Kamensky void *predictor_ctx; 44818efac9SLisandro Dalcin } TS_Alpha; 45818efac9SLisandro Dalcin 46220f924aSDavid Kamensky /*@C 47220f924aSDavid Kamensky TSAlpha2SetPredictor - sets the callback for computing a predictor (i.e., initial guess 48220f924aSDavid Kamensky for the nonlinear solver). 49220f924aSDavid Kamensky 50220f924aSDavid Kamensky Input Parameters: 51220f924aSDavid Kamensky + ts - timestepping context 52220f924aSDavid Kamensky . predictor - callback to set the predictor in each step 53220f924aSDavid Kamensky - ctx - the application context, which may be set to `NULL` if not used 54220f924aSDavid Kamensky 55220f924aSDavid Kamensky Level: intermediate 56220f924aSDavid Kamensky 57220f924aSDavid Kamensky Notes: 58220f924aSDavid Kamensky 59220f924aSDavid Kamensky If this function is never called, a same-state-vector predictor will be used, i.e., 60220f924aSDavid Kamensky the initial guess will be the converged solution from the previous time step, without regard 61220f924aSDavid Kamensky for the previous velocity or acceleration. 62220f924aSDavid Kamensky 638434afd1SBarry Smith .seealso: [](ch_ts), `TS`, `TSALPHA2`, `TSAlpha2PredictorFn` 64220f924aSDavid Kamensky @*/ 658434afd1SBarry Smith PetscErrorCode TSAlpha2SetPredictor(TS ts, TSAlpha2PredictorFn *predictor, void *ctx) 66220f924aSDavid Kamensky { 67*f4f49eeaSPierre Jolivet TS_Alpha *th = (TS_Alpha *)ts->data; 68220f924aSDavid Kamensky 69220f924aSDavid Kamensky PetscFunctionBegin; 70220f924aSDavid Kamensky th->predictor = predictor; 71220f924aSDavid Kamensky th->predictor_ctx = ctx; 72220f924aSDavid Kamensky PetscFunctionReturn(PETSC_SUCCESS); 73220f924aSDavid Kamensky } 74220f924aSDavid Kamensky 75220f924aSDavid Kamensky static PetscErrorCode TSAlpha_ApplyPredictor(TS ts, Vec X1) 76220f924aSDavid Kamensky { 77220f924aSDavid Kamensky /* Apply a custom predictor if set, or default to same-displacement. */ 78*f4f49eeaSPierre Jolivet TS_Alpha *th = (TS_Alpha *)ts->data; 79220f924aSDavid Kamensky 80220f924aSDavid Kamensky PetscFunctionBegin; 81220f924aSDavid Kamensky if (th->predictor) PetscCall(th->predictor(ts, th->X0, th->V0, th->A0, X1, th->predictor_ctx)); 82220f924aSDavid Kamensky else PetscCall(VecCopy(th->X0, X1)); 83220f924aSDavid Kamensky PetscFunctionReturn(PETSC_SUCCESS); 84220f924aSDavid Kamensky } 85220f924aSDavid 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)); 182220f924aSDavid 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)); 194220f924aSDavid 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: 222123b1737SDavid Kamensky /* Revert TSAlpha to the initial state (t0,X0,V0), but retain 223123b1737SDavid Kamensky potential time step reduction by factor after failed solve. */ 224818efac9SLisandro Dalcin if (initok) *initok = stageok; 225123b1737SDavid Kamensky PetscCall(TSSetTimeStep(ts, 2 * ts->time_step)); 2269566063dSJacob Faibussowitsch PetscCall(TSAlpha2SetParams(ts, alpha_m, alpha_f, gamma, beta)); 2279566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, th->X0)); 2289566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_dot, th->V0)); 229818efac9SLisandro Dalcin 2309566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X1)); 2319566063dSJacob Faibussowitsch PetscCall(VecDestroy(&V1)); 2323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 233818efac9SLisandro Dalcin } 234818efac9SLisandro Dalcin 235d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSStep_Alpha(TS ts) 236d71ae5a4SJacob Faibussowitsch { 237818efac9SLisandro Dalcin TS_Alpha *th = (TS_Alpha *)ts->data; 238818efac9SLisandro Dalcin PetscInt rejections = 0; 239818efac9SLisandro Dalcin PetscBool stageok, accept = PETSC_TRUE; 240818efac9SLisandro Dalcin PetscReal next_time_step = ts->time_step; 241818efac9SLisandro Dalcin 242818efac9SLisandro Dalcin PetscFunctionBegin; 2439566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister(citation, &cited)); 244818efac9SLisandro Dalcin 245818efac9SLisandro Dalcin if (!ts->steprollback) { 2469566063dSJacob Faibussowitsch if (th->vec_sol_prev) PetscCall(VecCopy(th->X0, th->vec_sol_prev)); 2479566063dSJacob Faibussowitsch if (th->vec_dot_prev) PetscCall(VecCopy(th->V0, th->vec_dot_prev)); 2489566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, th->X0)); 2499566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_dot, th->V0)); 2509566063dSJacob Faibussowitsch PetscCall(VecCopy(th->A1, th->A0)); 251818efac9SLisandro Dalcin } 252818efac9SLisandro Dalcin 253818efac9SLisandro Dalcin th->status = TS_STEP_INCOMPLETE; 254818efac9SLisandro Dalcin while (!ts->reason && th->status != TS_STEP_COMPLETE) { 255818efac9SLisandro Dalcin if (ts->steprestart) { 2569566063dSJacob Faibussowitsch PetscCall(TSAlpha_Restart(ts, &stageok)); 257818efac9SLisandro Dalcin if (!stageok) goto reject_step; 258818efac9SLisandro Dalcin } 259818efac9SLisandro Dalcin 2609566063dSJacob Faibussowitsch PetscCall(TSAlpha_StageTime(ts)); 261220f924aSDavid Kamensky PetscCall(TSAlpha_ApplyPredictor(ts, th->X1)); 2629566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts, th->stage_time)); 2639566063dSJacob Faibussowitsch PetscCall(TSAlpha_SNESSolve(ts, NULL, th->X1)); 2649566063dSJacob Faibussowitsch PetscCall(TSPostStage(ts, th->stage_time, 0, &th->Xa)); 2659566063dSJacob Faibussowitsch PetscCall(TSAdaptCheckStage(ts->adapt, ts, th->stage_time, th->Xa, &stageok)); 266818efac9SLisandro Dalcin if (!stageok) goto reject_step; 267818efac9SLisandro Dalcin 268818efac9SLisandro Dalcin th->status = TS_STEP_PENDING; 2699566063dSJacob Faibussowitsch PetscCall(VecCopy(th->X1, ts->vec_sol)); 2709566063dSJacob Faibussowitsch PetscCall(VecCopy(th->V1, ts->vec_dot)); 2719566063dSJacob Faibussowitsch PetscCall(TSAdaptChoose(ts->adapt, ts, ts->time_step, NULL, &next_time_step, &accept)); 272818efac9SLisandro Dalcin th->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 273818efac9SLisandro Dalcin if (!accept) { 2749566063dSJacob Faibussowitsch PetscCall(VecCopy(th->X0, ts->vec_sol)); 2759566063dSJacob Faibussowitsch PetscCall(VecCopy(th->V0, ts->vec_dot)); 276818efac9SLisandro Dalcin ts->time_step = next_time_step; 277818efac9SLisandro Dalcin goto reject_step; 278818efac9SLisandro Dalcin } 279818efac9SLisandro Dalcin 280818efac9SLisandro Dalcin ts->ptime += ts->time_step; 281818efac9SLisandro Dalcin ts->time_step = next_time_step; 282818efac9SLisandro Dalcin break; 283818efac9SLisandro Dalcin 284818efac9SLisandro Dalcin reject_step: 2859371c9d4SSatish Balay ts->reject++; 2869371c9d4SSatish Balay accept = PETSC_FALSE; 287818efac9SLisandro Dalcin if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 288818efac9SLisandro Dalcin ts->reason = TS_DIVERGED_STEP_REJECTED; 28963a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, rejections)); 290818efac9SLisandro Dalcin } 291818efac9SLisandro Dalcin } 2923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 293818efac9SLisandro Dalcin } 294818efac9SLisandro Dalcin 295d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSEvaluateWLTE_Alpha(TS ts, NormType wnormtype, PetscInt *order, PetscReal *wlte) 296d71ae5a4SJacob Faibussowitsch { 297818efac9SLisandro Dalcin TS_Alpha *th = (TS_Alpha *)ts->data; 298818efac9SLisandro Dalcin Vec X = th->X1; /* X = solution */ 299818efac9SLisandro Dalcin Vec V = th->V1; /* V = solution */ 300818efac9SLisandro Dalcin Vec Y = th->vec_lte_work[0]; /* Y = X + LTE */ 301818efac9SLisandro Dalcin Vec Z = th->vec_lte_work[1]; /* Z = V + LTE */ 3027453f775SEmil Constantinescu PetscReal enormX, enormV, enormXa, enormVa, enormXr, enormVr; 303818efac9SLisandro Dalcin 304818efac9SLisandro Dalcin PetscFunctionBegin; 3059371c9d4SSatish Balay if (!th->vec_sol_prev) { 3069371c9d4SSatish Balay *wlte = -1; 3073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3089371c9d4SSatish Balay } 3099371c9d4SSatish Balay if (!th->vec_dot_prev) { 3109371c9d4SSatish Balay *wlte = -1; 3113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3129371c9d4SSatish Balay } 3139371c9d4SSatish Balay if (!th->vec_lte_work[0]) { 3149371c9d4SSatish Balay *wlte = -1; 3153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3169371c9d4SSatish Balay } 3179371c9d4SSatish Balay if (!th->vec_lte_work[1]) { 3189371c9d4SSatish Balay *wlte = -1; 3193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3209371c9d4SSatish Balay } 321818efac9SLisandro Dalcin if (ts->steprestart) { 3222ffb9264SLisandro Dalcin /* th->vec_lte_prev is set to the LTE in TSAlpha_Restart() */ 3239566063dSJacob Faibussowitsch PetscCall(VecAXPY(Y, 1, X)); 3249566063dSJacob Faibussowitsch PetscCall(VecAXPY(Z, 1, V)); 325818efac9SLisandro Dalcin } else { 326818efac9SLisandro Dalcin /* Compute LTE using backward differences with non-constant time step */ 327818efac9SLisandro Dalcin PetscReal h = ts->time_step, h_prev = ts->ptime - ts->ptime_prev; 328818efac9SLisandro Dalcin PetscReal a = 1 + h_prev / h; 3299371c9d4SSatish Balay PetscScalar scal[3]; 3309371c9d4SSatish Balay Vec vecX[3], vecV[3]; 3319371c9d4SSatish Balay scal[0] = +1 / a; 3329371c9d4SSatish Balay scal[1] = -1 / (a - 1); 3339371c9d4SSatish Balay scal[2] = +1 / (a * (a - 1)); 3349371c9d4SSatish Balay vecX[0] = th->X1; 3359371c9d4SSatish Balay vecX[1] = th->X0; 3369371c9d4SSatish Balay vecX[2] = th->vec_sol_prev; 3379371c9d4SSatish Balay vecV[0] = th->V1; 3389371c9d4SSatish Balay vecV[1] = th->V0; 3399371c9d4SSatish Balay vecV[2] = th->vec_dot_prev; 3409566063dSJacob Faibussowitsch PetscCall(VecCopy(X, Y)); 3419566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Y, 3, scal, vecX)); 3429566063dSJacob Faibussowitsch PetscCall(VecCopy(V, Z)); 3439566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Z, 3, scal, vecV)); 344818efac9SLisandro Dalcin } 345818efac9SLisandro Dalcin /* XXX ts->atol and ts->vatol are not appropriate for computing enormV */ 3469566063dSJacob Faibussowitsch PetscCall(TSErrorWeightedNorm(ts, X, Y, wnormtype, &enormX, &enormXa, &enormXr)); 3479566063dSJacob Faibussowitsch PetscCall(TSErrorWeightedNorm(ts, V, Z, wnormtype, &enormV, &enormVa, &enormVr)); 3489371c9d4SSatish Balay if (wnormtype == NORM_2) *wlte = PetscSqrtReal(PetscSqr(enormX) / 2 + PetscSqr(enormV) / 2); 3499371c9d4SSatish Balay else *wlte = PetscMax(enormX, enormV); 350818efac9SLisandro Dalcin if (order) *order = 2; 3513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 352818efac9SLisandro Dalcin } 353818efac9SLisandro Dalcin 354d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRollBack_Alpha(TS ts) 355d71ae5a4SJacob Faibussowitsch { 356818efac9SLisandro Dalcin TS_Alpha *th = (TS_Alpha *)ts->data; 357818efac9SLisandro Dalcin 358818efac9SLisandro Dalcin PetscFunctionBegin; 3599566063dSJacob Faibussowitsch PetscCall(VecCopy(th->X0, ts->vec_sol)); 3609566063dSJacob Faibussowitsch PetscCall(VecCopy(th->V0, ts->vec_dot)); 3613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 362818efac9SLisandro Dalcin } 363818efac9SLisandro Dalcin 364818efac9SLisandro Dalcin /* 365818efac9SLisandro Dalcin static PetscErrorCode TSInterpolate_Alpha(TS ts,PetscReal t,Vec X,Vec V) 366818efac9SLisandro Dalcin { 367818efac9SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 368818efac9SLisandro Dalcin PetscReal dt = t - ts->ptime; 369818efac9SLisandro Dalcin 370818efac9SLisandro Dalcin PetscFunctionBegin; 3719566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_dot,V)); 3729566063dSJacob Faibussowitsch PetscCall(VecAXPY(V,dt*(1-th->Gamma),th->A0)); 3739566063dSJacob Faibussowitsch PetscCall(VecAXPY(V,dt*th->Gamma,th->A1)); 3749566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol,X)); 3759566063dSJacob Faibussowitsch PetscCall(VecAXPY(X,dt,V)); 3769566063dSJacob Faibussowitsch PetscCall(VecAXPY(X,dt*dt*((PetscReal)0.5-th->Beta),th->A0)); 3779566063dSJacob Faibussowitsch PetscCall(VecAXPY(X,dt*dt*th->Beta,th->A1)); 3783ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 379818efac9SLisandro Dalcin } 380818efac9SLisandro Dalcin */ 381818efac9SLisandro Dalcin 382d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormFunction_Alpha(PETSC_UNUSED SNES snes, Vec X, Vec F, TS ts) 383d71ae5a4SJacob Faibussowitsch { 384818efac9SLisandro Dalcin TS_Alpha *th = (TS_Alpha *)ts->data; 385818efac9SLisandro Dalcin PetscReal ta = th->stage_time; 386818efac9SLisandro Dalcin Vec Xa = th->Xa, Va = th->Va, Aa = th->Aa; 387818efac9SLisandro Dalcin 388818efac9SLisandro Dalcin PetscFunctionBegin; 3899566063dSJacob Faibussowitsch PetscCall(TSAlpha_StageVecs(ts, X)); 390818efac9SLisandro Dalcin /* F = Function(ta,Xa,Va,Aa) */ 3919566063dSJacob Faibussowitsch PetscCall(TSComputeI2Function(ts, ta, Xa, Va, Aa, F)); 3929566063dSJacob Faibussowitsch PetscCall(VecScale(F, th->scale_F)); 3933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 394818efac9SLisandro Dalcin } 395818efac9SLisandro Dalcin 396d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormJacobian_Alpha(PETSC_UNUSED SNES snes, PETSC_UNUSED Vec X, Mat J, Mat P, TS ts) 397d71ae5a4SJacob Faibussowitsch { 398818efac9SLisandro Dalcin TS_Alpha *th = (TS_Alpha *)ts->data; 399818efac9SLisandro Dalcin PetscReal ta = th->stage_time; 400818efac9SLisandro Dalcin Vec Xa = th->Xa, Va = th->Va, Aa = th->Aa; 401818efac9SLisandro Dalcin PetscReal dVdX = th->shift_V, dAdX = th->shift_A; 402818efac9SLisandro Dalcin 403818efac9SLisandro Dalcin PetscFunctionBegin; 404818efac9SLisandro Dalcin /* J,P = Jacobian(ta,Xa,Va,Aa) */ 4059566063dSJacob Faibussowitsch PetscCall(TSComputeI2Jacobian(ts, ta, Xa, Va, Aa, dVdX, dAdX, J, P)); 4063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 407818efac9SLisandro Dalcin } 408818efac9SLisandro Dalcin 409d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSReset_Alpha(TS ts) 410d71ae5a4SJacob Faibussowitsch { 411818efac9SLisandro Dalcin TS_Alpha *th = (TS_Alpha *)ts->data; 412818efac9SLisandro Dalcin 413818efac9SLisandro Dalcin PetscFunctionBegin; 4149566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->X0)); 4159566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->Xa)); 4169566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->X1)); 4179566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->V0)); 4189566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->Va)); 4199566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->V1)); 4209566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->A0)); 4219566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->Aa)); 4229566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->A1)); 4239566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->vec_sol_prev)); 4249566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->vec_dot_prev)); 4259566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->vec_lte_work[0])); 4269566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->vec_lte_work[1])); 4273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 428818efac9SLisandro Dalcin } 429818efac9SLisandro Dalcin 430d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSDestroy_Alpha(TS ts) 431d71ae5a4SJacob Faibussowitsch { 432818efac9SLisandro Dalcin PetscFunctionBegin; 4339566063dSJacob Faibussowitsch PetscCall(TSReset_Alpha(ts)); 4349566063dSJacob Faibussowitsch PetscCall(PetscFree(ts->data)); 435818efac9SLisandro Dalcin 4369566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlpha2SetRadius_C", NULL)); 4379566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlpha2SetParams_C", NULL)); 4389566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlpha2GetParams_C", NULL)); 4393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 440818efac9SLisandro Dalcin } 441818efac9SLisandro Dalcin 442d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_Alpha(TS ts) 443d71ae5a4SJacob Faibussowitsch { 444818efac9SLisandro Dalcin TS_Alpha *th = (TS_Alpha *)ts->data; 4452ffb9264SLisandro Dalcin PetscBool match; 446818efac9SLisandro Dalcin 447818efac9SLisandro Dalcin PetscFunctionBegin; 4489566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &th->X0)); 4499566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &th->Xa)); 4509566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &th->X1)); 4519566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &th->V0)); 4529566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &th->Va)); 4539566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &th->V1)); 4549566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &th->A0)); 4559566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &th->Aa)); 4569566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &th->A1)); 457818efac9SLisandro Dalcin 4589566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts, &ts->adapt)); 4599566063dSJacob Faibussowitsch PetscCall(TSAdaptCandidatesClear(ts->adapt)); 4609566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)ts->adapt, TSADAPTNONE, &match)); 4612ffb9264SLisandro Dalcin if (!match) { 4629566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &th->vec_sol_prev)); 4639566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &th->vec_dot_prev)); 4649566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &th->vec_lte_work[0])); 4659566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &th->vec_lte_work[1])); 466818efac9SLisandro Dalcin } 467818efac9SLisandro Dalcin 4689566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &ts->snes)); 4693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 470818efac9SLisandro Dalcin } 471818efac9SLisandro Dalcin 472d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetFromOptions_Alpha(TS ts, PetscOptionItems *PetscOptionsObject) 473d71ae5a4SJacob Faibussowitsch { 474818efac9SLisandro Dalcin TS_Alpha *th = (TS_Alpha *)ts->data; 475818efac9SLisandro Dalcin 476818efac9SLisandro Dalcin PetscFunctionBegin; 477d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Generalized-Alpha ODE solver options"); 478818efac9SLisandro Dalcin { 479818efac9SLisandro Dalcin PetscBool flg; 480818efac9SLisandro Dalcin PetscReal radius = 1; 4819566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_alpha_radius", "Spectral radius (high-frequency dissipation)", "TSAlpha2SetRadius", radius, &radius, &flg)); 4829566063dSJacob Faibussowitsch if (flg) PetscCall(TSAlpha2SetRadius(ts, radius)); 4839566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_alpha_alpha_m", "Algorithmic parameter alpha_m", "TSAlpha2SetParams", th->Alpha_m, &th->Alpha_m, NULL)); 4849566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_alpha_alpha_f", "Algorithmic parameter alpha_f", "TSAlpha2SetParams", th->Alpha_f, &th->Alpha_f, NULL)); 4859566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_alpha_gamma", "Algorithmic parameter gamma", "TSAlpha2SetParams", th->Gamma, &th->Gamma, NULL)); 4869566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_alpha_beta", "Algorithmic parameter beta", "TSAlpha2SetParams", th->Beta, &th->Beta, NULL)); 4879566063dSJacob Faibussowitsch PetscCall(TSAlpha2SetParams(ts, th->Alpha_m, th->Alpha_f, th->Gamma, th->Beta)); 488818efac9SLisandro Dalcin } 489d0609cedSBarry Smith PetscOptionsHeadEnd(); 4903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 491818efac9SLisandro Dalcin } 492818efac9SLisandro Dalcin 493d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSView_Alpha(TS ts, PetscViewer viewer) 494d71ae5a4SJacob Faibussowitsch { 495818efac9SLisandro Dalcin TS_Alpha *th = (TS_Alpha *)ts->data; 496818efac9SLisandro Dalcin PetscBool iascii; 497818efac9SLisandro Dalcin 498818efac9SLisandro Dalcin PetscFunctionBegin; 4999566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 50048a46eb9SPierre 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)); 5013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 502818efac9SLisandro Dalcin } 503818efac9SLisandro Dalcin 504d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAlpha2SetRadius_Alpha(TS ts, PetscReal radius) 505d71ae5a4SJacob Faibussowitsch { 506818efac9SLisandro Dalcin PetscReal alpha_m, alpha_f, gamma, beta; 507818efac9SLisandro Dalcin 508818efac9SLisandro Dalcin PetscFunctionBegin; 509cad9d221SBarry Smith PetscCheck(radius >= 0 && radius <= 1, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_OUTOFRANGE, "Radius %g not in range [0,1]", (double)radius); 510818efac9SLisandro Dalcin alpha_m = (2 - radius) / (1 + radius); 511818efac9SLisandro Dalcin alpha_f = 1 / (1 + radius); 512818efac9SLisandro Dalcin gamma = (PetscReal)0.5 + alpha_m - alpha_f; 5139371c9d4SSatish Balay beta = (PetscReal)0.5 * (1 + alpha_m - alpha_f); 5149371c9d4SSatish Balay beta *= beta; 5159566063dSJacob Faibussowitsch PetscCall(TSAlpha2SetParams(ts, alpha_m, alpha_f, gamma, beta)); 5163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 517818efac9SLisandro Dalcin } 518818efac9SLisandro Dalcin 519d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAlpha2SetParams_Alpha(TS ts, PetscReal alpha_m, PetscReal alpha_f, PetscReal gamma, PetscReal beta) 520d71ae5a4SJacob Faibussowitsch { 521818efac9SLisandro Dalcin TS_Alpha *th = (TS_Alpha *)ts->data; 522818efac9SLisandro Dalcin PetscReal tol = 100 * PETSC_MACHINE_EPSILON; 523818efac9SLisandro Dalcin PetscReal res = ((PetscReal)0.5 + alpha_m - alpha_f) - gamma; 524818efac9SLisandro Dalcin 525818efac9SLisandro Dalcin PetscFunctionBegin; 526818efac9SLisandro Dalcin th->Alpha_m = alpha_m; 527818efac9SLisandro Dalcin th->Alpha_f = alpha_f; 528818efac9SLisandro Dalcin th->Gamma = gamma; 529818efac9SLisandro Dalcin th->Beta = beta; 530818efac9SLisandro Dalcin th->order = (PetscAbsReal(res) < tol) ? 2 : 1; 5313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 532818efac9SLisandro Dalcin } 533818efac9SLisandro Dalcin 534d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAlpha2GetParams_Alpha(TS ts, PetscReal *alpha_m, PetscReal *alpha_f, PetscReal *gamma, PetscReal *beta) 535d71ae5a4SJacob Faibussowitsch { 536818efac9SLisandro Dalcin TS_Alpha *th = (TS_Alpha *)ts->data; 537818efac9SLisandro Dalcin 538818efac9SLisandro Dalcin PetscFunctionBegin; 539818efac9SLisandro Dalcin if (alpha_m) *alpha_m = th->Alpha_m; 540818efac9SLisandro Dalcin if (alpha_f) *alpha_f = th->Alpha_f; 541818efac9SLisandro Dalcin if (gamma) *gamma = th->Gamma; 542818efac9SLisandro Dalcin if (beta) *beta = th->Beta; 5433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 544818efac9SLisandro Dalcin } 545818efac9SLisandro Dalcin 546818efac9SLisandro Dalcin /*MC 5471d27aa22SBarry Smith TSALPHA2 - ODE/DAE solver using the implicit Generalized-Alpha method for second-order systems {cite}`chung1993` 548818efac9SLisandro Dalcin 549818efac9SLisandro Dalcin Level: beginner 550818efac9SLisandro Dalcin 5511cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSCreate()`, `TSSetType()`, `TSAlpha2SetRadius()`, `TSAlpha2SetParams()` 552818efac9SLisandro Dalcin M*/ 553d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_Alpha2(TS ts) 554d71ae5a4SJacob Faibussowitsch { 555818efac9SLisandro Dalcin TS_Alpha *th; 556818efac9SLisandro Dalcin 557818efac9SLisandro Dalcin PetscFunctionBegin; 558818efac9SLisandro Dalcin ts->ops->reset = TSReset_Alpha; 559818efac9SLisandro Dalcin ts->ops->destroy = TSDestroy_Alpha; 560818efac9SLisandro Dalcin ts->ops->view = TSView_Alpha; 561818efac9SLisandro Dalcin ts->ops->setup = TSSetUp_Alpha; 562818efac9SLisandro Dalcin ts->ops->setfromoptions = TSSetFromOptions_Alpha; 563818efac9SLisandro Dalcin ts->ops->step = TSStep_Alpha; 564818efac9SLisandro Dalcin ts->ops->evaluatewlte = TSEvaluateWLTE_Alpha; 565818efac9SLisandro Dalcin ts->ops->rollback = TSRollBack_Alpha; 566818efac9SLisandro Dalcin /*ts->ops->interpolate = TSInterpolate_Alpha;*/ 567818efac9SLisandro Dalcin ts->ops->snesfunction = SNESTSFormFunction_Alpha; 568818efac9SLisandro Dalcin ts->ops->snesjacobian = SNESTSFormJacobian_Alpha; 5692ffb9264SLisandro Dalcin ts->default_adapt_type = TSADAPTNONE; 570818efac9SLisandro Dalcin 571efd4aadfSBarry Smith ts->usessnes = PETSC_TRUE; 572efd4aadfSBarry Smith 5734dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&th)); 574818efac9SLisandro Dalcin ts->data = (void *)th; 575818efac9SLisandro Dalcin 576818efac9SLisandro Dalcin th->Alpha_m = 0.5; 577818efac9SLisandro Dalcin th->Alpha_f = 0.5; 578818efac9SLisandro Dalcin th->Gamma = 0.5; 579818efac9SLisandro Dalcin th->Beta = 0.25; 580818efac9SLisandro Dalcin th->order = 2; 581818efac9SLisandro Dalcin 582220f924aSDavid Kamensky th->predictor = NULL; 583220f924aSDavid Kamensky th->predictor_ctx = NULL; 584220f924aSDavid Kamensky 5859566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlpha2SetRadius_C", TSAlpha2SetRadius_Alpha)); 5869566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlpha2SetParams_C", TSAlpha2SetParams_Alpha)); 5879566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlpha2GetParams_C", TSAlpha2GetParams_Alpha)); 5883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 589818efac9SLisandro Dalcin } 590818efac9SLisandro Dalcin 591818efac9SLisandro Dalcin /*@ 592bcf0153eSBarry Smith TSAlpha2SetRadius - sets the desired spectral radius of the method for `TSALPHA2` 593818efac9SLisandro Dalcin (i.e. high-frequency numerical damping) 594818efac9SLisandro Dalcin 595c3339decSBarry Smith Logically Collective 596818efac9SLisandro Dalcin 597d8d19677SJose E. Roman Input Parameters: 598818efac9SLisandro Dalcin + ts - timestepping context 599818efac9SLisandro Dalcin - radius - the desired spectral radius 600818efac9SLisandro Dalcin 601bcf0153eSBarry Smith Options Database Key: 60267b8a455SSatish Balay . -ts_alpha_radius <radius> - set the desired spectral radius 603818efac9SLisandro Dalcin 604818efac9SLisandro Dalcin Level: intermediate 605818efac9SLisandro Dalcin 60614d0ab18SJacob Faibussowitsch Notes: 60714d0ab18SJacob Faibussowitsch 60814d0ab18SJacob Faibussowitsch The algorithmic parameters $\alpha_m$ and $\alpha_f$ of the generalized-$\alpha$ method can 60914d0ab18SJacob Faibussowitsch be computed in terms of a specified spectral radius $\rho$ in `[0, 1]` for infinite time step 61014d0ab18SJacob Faibussowitsch in order to control high-frequency numerical damping\: 611562efe2eSBarry Smith 61214d0ab18SJacob Faibussowitsch $$ 613562efe2eSBarry Smith \begin{align*} 614562efe2eSBarry Smith \alpha_m = (2-\rho)/(1+\rho) \\ 61514d0ab18SJacob Faibussowitsch \alpha_f = 1/(1+\rho) 616562efe2eSBarry Smith \end{align*} 61714d0ab18SJacob Faibussowitsch $$ 61814d0ab18SJacob Faibussowitsch 6191cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSALPHA2`, `TSAlpha2SetParams()`, `TSAlpha2GetParams()` 620818efac9SLisandro Dalcin @*/ 621d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAlpha2SetRadius(TS ts, PetscReal radius) 622d71ae5a4SJacob Faibussowitsch { 623818efac9SLisandro Dalcin PetscFunctionBegin; 624818efac9SLisandro Dalcin PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 625818efac9SLisandro Dalcin PetscValidLogicalCollectiveReal(ts, radius, 2); 626cad9d221SBarry Smith PetscCheck(radius >= 0 && radius <= 1, ((PetscObject)ts)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Radius %g not in range [0,1]", (double)radius); 627cac4c232SBarry Smith PetscTryMethod(ts, "TSAlpha2SetRadius_C", (TS, PetscReal), (ts, radius)); 6283ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 629818efac9SLisandro Dalcin } 630818efac9SLisandro Dalcin 631818efac9SLisandro Dalcin /*@ 632bcf0153eSBarry Smith TSAlpha2SetParams - sets the algorithmic parameters for `TSALPHA2` 633818efac9SLisandro Dalcin 634c3339decSBarry Smith Logically Collective 635818efac9SLisandro Dalcin 636d8d19677SJose E. Roman Input Parameters: 637818efac9SLisandro Dalcin + ts - timestepping context 6382fe279fdSBarry Smith . alpha_m - algorithmic parameter 6392fe279fdSBarry Smith . alpha_f - algorithmic parameter 6402fe279fdSBarry Smith . gamma - algorithmic parameter 6412fe279fdSBarry Smith - beta - algorithmic parameter 642818efac9SLisandro Dalcin 643bcf0153eSBarry Smith Options Database Keys: 64467b8a455SSatish Balay + -ts_alpha_alpha_m <alpha_m> - set alpha_m 64567b8a455SSatish Balay . -ts_alpha_alpha_f <alpha_f> - set alpha_f 64667b8a455SSatish Balay . -ts_alpha_gamma <gamma> - set gamma 64767b8a455SSatish Balay - -ts_alpha_beta <beta> - set beta 648818efac9SLisandro Dalcin 649bcf0153eSBarry Smith Level: advanced 650bcf0153eSBarry Smith 65114d0ab18SJacob Faibussowitsch Notes: 65214d0ab18SJacob Faibussowitsch Second-order accuracy can be obtained so long as\: 653562efe2eSBarry Smith 65414d0ab18SJacob Faibussowitsch $$ 655562efe2eSBarry Smith \begin{align*} 656562efe2eSBarry Smith \gamma = 1/2 + \alpha_m - \alpha_f \\ 657562efe2eSBarry Smith \beta = 1/4 (1 + \alpha_m - \alpha_f)^2. 658562efe2eSBarry Smith \end{align*} 65914d0ab18SJacob Faibussowitsch $$ 66014d0ab18SJacob Faibussowitsch 66114d0ab18SJacob Faibussowitsch Unconditional stability requires\: 66214d0ab18SJacob Faibussowitsch $$ 663562efe2eSBarry Smith \alpha_m >= \alpha_f >= 1/2. 66414d0ab18SJacob Faibussowitsch $$ 66514d0ab18SJacob Faibussowitsch 66614d0ab18SJacob Faibussowitsch Use of this function is normally only required to hack `TSALPHA2` to use a modified 66714d0ab18SJacob Faibussowitsch integration scheme. Users should call `TSAlpha2SetRadius()` to set the desired spectral 66814d0ab18SJacob Faibussowitsch radius of the methods (i.e. high-frequency damping) in order so select optimal values for 669818efac9SLisandro Dalcin these parameters. 670818efac9SLisandro Dalcin 6711cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSALPHA2`, `TSAlpha2SetRadius()`, `TSAlpha2GetParams()` 672818efac9SLisandro Dalcin @*/ 673d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAlpha2SetParams(TS ts, PetscReal alpha_m, PetscReal alpha_f, PetscReal gamma, PetscReal beta) 674d71ae5a4SJacob Faibussowitsch { 675818efac9SLisandro Dalcin PetscFunctionBegin; 676818efac9SLisandro Dalcin PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 677818efac9SLisandro Dalcin PetscValidLogicalCollectiveReal(ts, alpha_m, 2); 678818efac9SLisandro Dalcin PetscValidLogicalCollectiveReal(ts, alpha_f, 3); 679818efac9SLisandro Dalcin PetscValidLogicalCollectiveReal(ts, gamma, 4); 680818efac9SLisandro Dalcin PetscValidLogicalCollectiveReal(ts, beta, 5); 681cac4c232SBarry Smith PetscTryMethod(ts, "TSAlpha2SetParams_C", (TS, PetscReal, PetscReal, PetscReal, PetscReal), (ts, alpha_m, alpha_f, gamma, beta)); 6823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 683818efac9SLisandro Dalcin } 684818efac9SLisandro Dalcin 685818efac9SLisandro Dalcin /*@ 686bcf0153eSBarry Smith TSAlpha2GetParams - gets the algorithmic parameters for `TSALPHA2` 687818efac9SLisandro Dalcin 688818efac9SLisandro Dalcin Not Collective 689818efac9SLisandro Dalcin 690818efac9SLisandro Dalcin Input Parameter: 691818efac9SLisandro Dalcin . ts - timestepping context 692818efac9SLisandro Dalcin 693818efac9SLisandro Dalcin Output Parameters: 6942fe279fdSBarry Smith + alpha_m - algorithmic parameter 6952fe279fdSBarry Smith . alpha_f - algorithmic parameter 6962fe279fdSBarry Smith . gamma - algorithmic parameter 6972fe279fdSBarry Smith - beta - algorithmic parameter 698818efac9SLisandro Dalcin 699bcf0153eSBarry Smith Level: advanced 700bcf0153eSBarry Smith 701818efac9SLisandro Dalcin Note: 70214d0ab18SJacob Faibussowitsch Use of this function is normally only required to hack `TSALPHA2` to use a modified 70314d0ab18SJacob Faibussowitsch integration scheme. Users should call `TSAlpha2SetRadius()` to set the high-frequency damping 70414d0ab18SJacob Faibussowitsch (i.e. spectral radius of the method) in order so select optimal values for these parameters. 705818efac9SLisandro Dalcin 7061cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSALPHA2`, `TSAlpha2SetRadius()`, `TSAlpha2SetParams()` 707818efac9SLisandro Dalcin @*/ 708d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAlpha2GetParams(TS ts, PetscReal *alpha_m, PetscReal *alpha_f, PetscReal *gamma, PetscReal *beta) 709d71ae5a4SJacob Faibussowitsch { 710818efac9SLisandro Dalcin PetscFunctionBegin; 711818efac9SLisandro Dalcin PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 7124f572ea9SToby Isaac if (alpha_m) PetscAssertPointer(alpha_m, 2); 7134f572ea9SToby Isaac if (alpha_f) PetscAssertPointer(alpha_f, 3); 7144f572ea9SToby Isaac if (gamma) PetscAssertPointer(gamma, 4); 7154f572ea9SToby Isaac if (beta) PetscAssertPointer(beta, 5); 716cac4c232SBarry Smith PetscUseMethod(ts, "TSAlpha2GetParams_C", (TS, PetscReal *, PetscReal *, PetscReal *, PetscReal *), (ts, alpha_m, alpha_f, gamma, beta)); 7173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 718818efac9SLisandro Dalcin } 719