1b07a2398SLisandro Dalcin /* 2b07a2398SLisandro Dalcin Code for timestepping with implicit generalized-\alpha method 3b07a2398SLisandro Dalcin for first order systems. 4b07a2398SLisandro Dalcin */ 5b07a2398SLisandro Dalcin #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 6b07a2398SLisandro Dalcin 7b07a2398SLisandro Dalcin static PetscBool cited = PETSC_FALSE; 8b07a2398SLisandro Dalcin static const char citation[] = 9b07a2398SLisandro Dalcin "@article{Jansen2000,\n" 10b07a2398SLisandro Dalcin " title = {A generalized-$\\alpha$ method for integrating the filtered {N}avier--{S}tokes equations with a stabilized finite element method},\n" 11b07a2398SLisandro Dalcin " author = {Kenneth E. Jansen and Christian H. Whiting and Gregory M. Hulbert},\n" 12b07a2398SLisandro Dalcin " journal = {Computer Methods in Applied Mechanics and Engineering},\n" 13b07a2398SLisandro Dalcin " volume = {190},\n" 14b07a2398SLisandro Dalcin " number = {3--4},\n" 15b07a2398SLisandro Dalcin " pages = {305--319},\n" 16b07a2398SLisandro Dalcin " year = {2000},\n" 17b07a2398SLisandro Dalcin " issn = {0045-7825},\n" 18b07a2398SLisandro Dalcin " doi = {http://dx.doi.org/10.1016/S0045-7825(00)00203-6}\n}\n"; 19b07a2398SLisandro Dalcin 20b07a2398SLisandro Dalcin typedef struct { 21b07a2398SLisandro Dalcin PetscReal stage_time; 22b07a2398SLisandro Dalcin PetscReal shift_V; 23b07a2398SLisandro Dalcin PetscReal scale_F; 24b07a2398SLisandro Dalcin Vec X0,Xa,X1; 25b07a2398SLisandro Dalcin Vec V0,Va,V1; 26b07a2398SLisandro Dalcin 27b07a2398SLisandro Dalcin PetscReal Alpha_m; 28b07a2398SLisandro Dalcin PetscReal Alpha_f; 29b07a2398SLisandro Dalcin PetscReal Gamma; 30b07a2398SLisandro Dalcin PetscInt order; 31b07a2398SLisandro Dalcin 32b07a2398SLisandro Dalcin Vec vec_sol_prev; 331566a47fSLisandro Dalcin Vec vec_lte_work; 34b07a2398SLisandro Dalcin 35b07a2398SLisandro Dalcin TSStepStatus status; 36b07a2398SLisandro Dalcin } TS_Alpha; 37b07a2398SLisandro Dalcin 38b07a2398SLisandro Dalcin static PetscErrorCode TSAlpha_StageTime(TS ts) 39b07a2398SLisandro Dalcin { 40b07a2398SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 41b07a2398SLisandro Dalcin PetscReal t = ts->ptime; 42b07a2398SLisandro Dalcin PetscReal dt = ts->time_step; 43b07a2398SLisandro Dalcin PetscReal Alpha_m = th->Alpha_m; 44b07a2398SLisandro Dalcin PetscReal Alpha_f = th->Alpha_f; 45b07a2398SLisandro Dalcin PetscReal Gamma = th->Gamma; 46b07a2398SLisandro Dalcin 47b07a2398SLisandro Dalcin PetscFunctionBegin; 48b07a2398SLisandro Dalcin th->stage_time = t + Alpha_f*dt; 49b07a2398SLisandro Dalcin th->shift_V = Alpha_m/(Alpha_f*Gamma*dt); 50b07a2398SLisandro Dalcin th->scale_F = 1/Alpha_f; 51b07a2398SLisandro Dalcin PetscFunctionReturn(0); 52b07a2398SLisandro Dalcin } 53b07a2398SLisandro Dalcin 54b07a2398SLisandro Dalcin static PetscErrorCode TSAlpha_StageVecs(TS ts,Vec X) 55b07a2398SLisandro Dalcin { 56b07a2398SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 57b07a2398SLisandro Dalcin Vec X1 = X, V1 = th->V1; 58b07a2398SLisandro Dalcin Vec Xa = th->Xa, Va = th->Va; 59b07a2398SLisandro Dalcin Vec X0 = th->X0, V0 = th->V0; 60b07a2398SLisandro Dalcin PetscReal dt = ts->time_step; 61b07a2398SLisandro Dalcin PetscReal Alpha_m = th->Alpha_m; 62b07a2398SLisandro Dalcin PetscReal Alpha_f = th->Alpha_f; 63b07a2398SLisandro Dalcin PetscReal Gamma = th->Gamma; 64b07a2398SLisandro Dalcin 65b07a2398SLisandro Dalcin PetscFunctionBegin; 66b07a2398SLisandro Dalcin /* V1 = 1/(Gamma*dT)*(X1-X0) + (1-1/Gamma)*V0 */ 679566063dSJacob Faibussowitsch PetscCall(VecWAXPY(V1,-1.0,X0,X1)); 689566063dSJacob Faibussowitsch PetscCall(VecAXPBY(V1,1-1/Gamma,1/(Gamma*dt),V0)); 69b07a2398SLisandro Dalcin /* Xa = X0 + Alpha_f*(X1-X0) */ 709566063dSJacob Faibussowitsch PetscCall(VecWAXPY(Xa,-1.0,X0,X1)); 719566063dSJacob Faibussowitsch PetscCall(VecAYPX(Xa,Alpha_f,X0)); 72b07a2398SLisandro Dalcin /* Va = V0 + Alpha_m*(V1-V0) */ 739566063dSJacob Faibussowitsch PetscCall(VecWAXPY(Va,-1.0,V0,V1)); 749566063dSJacob Faibussowitsch PetscCall(VecAYPX(Va,Alpha_m,V0)); 75b07a2398SLisandro Dalcin PetscFunctionReturn(0); 76b07a2398SLisandro Dalcin } 77b07a2398SLisandro Dalcin 78470880abSPatrick Sanan static PetscErrorCode TSAlpha_SNESSolve(TS ts,Vec b,Vec x) 79b07a2398SLisandro Dalcin { 80b07a2398SLisandro Dalcin PetscInt nits,lits; 81b07a2398SLisandro Dalcin 82b07a2398SLisandro Dalcin PetscFunctionBegin; 839566063dSJacob Faibussowitsch PetscCall(SNESSolve(ts->snes,b,x)); 849566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(ts->snes,&nits)); 859566063dSJacob Faibussowitsch PetscCall(SNESGetLinearSolveIterations(ts->snes,&lits)); 86b07a2398SLisandro Dalcin ts->snes_its += nits; ts->ksp_its += lits; 87b07a2398SLisandro Dalcin PetscFunctionReturn(0); 88b07a2398SLisandro Dalcin } 89b07a2398SLisandro Dalcin 90b07a2398SLisandro Dalcin /* 91b07a2398SLisandro Dalcin Compute a consistent initial state for the generalized-alpha method. 92b07a2398SLisandro Dalcin - Solve two successive backward Euler steps with halved time step. 93b07a2398SLisandro Dalcin - Compute the initial time derivative using backward differences. 94b07a2398SLisandro Dalcin - If using adaptivity, estimate the LTE of the initial step. 95b07a2398SLisandro Dalcin */ 96fecfb714SLisandro Dalcin static PetscErrorCode TSAlpha_Restart(TS ts,PetscBool *initok) 97b07a2398SLisandro Dalcin { 98b07a2398SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 99b07a2398SLisandro Dalcin PetscReal time_step; 100b07a2398SLisandro Dalcin PetscReal alpha_m,alpha_f,gamma; 101b07a2398SLisandro Dalcin Vec X0 = ts->vec_sol, X1, X2 = th->X1; 102b07a2398SLisandro Dalcin PetscBool stageok; 103b07a2398SLisandro Dalcin 104b07a2398SLisandro Dalcin PetscFunctionBegin; 1059566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X0,&X1)); 106b07a2398SLisandro Dalcin 107b07a2398SLisandro Dalcin /* Setup backward Euler with halved time step */ 1089566063dSJacob Faibussowitsch PetscCall(TSAlphaGetParams(ts,&alpha_m,&alpha_f,&gamma)); 1099566063dSJacob Faibussowitsch PetscCall(TSAlphaSetParams(ts,1,1,1)); 1109566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts,&time_step)); 111b07a2398SLisandro Dalcin ts->time_step = time_step/2; 1129566063dSJacob Faibussowitsch PetscCall(TSAlpha_StageTime(ts)); 113b07a2398SLisandro Dalcin th->stage_time = ts->ptime; 1149566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(th->V0)); 115b07a2398SLisandro Dalcin 116b07a2398SLisandro Dalcin /* First BE step, (t0,X0) -> (t1,X1) */ 117b07a2398SLisandro Dalcin th->stage_time += ts->time_step; 1189566063dSJacob Faibussowitsch PetscCall(VecCopy(X0,th->X0)); 1199566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts,th->stage_time)); 1209566063dSJacob Faibussowitsch PetscCall(VecCopy(th->X0,X1)); 1219566063dSJacob Faibussowitsch PetscCall(TSAlpha_SNESSolve(ts,NULL,X1)); 1229566063dSJacob Faibussowitsch PetscCall(TSPostStage(ts,th->stage_time,0,&X1)); 1239566063dSJacob Faibussowitsch PetscCall(TSAdaptCheckStage(ts->adapt,ts,th->stage_time,X1,&stageok)); 124b07a2398SLisandro Dalcin if (!stageok) goto finally; 125b07a2398SLisandro Dalcin 126b07a2398SLisandro Dalcin /* Second BE step, (t1,X1) -> (t2,X2) */ 127b07a2398SLisandro Dalcin th->stage_time += ts->time_step; 1289566063dSJacob Faibussowitsch PetscCall(VecCopy(X1,th->X0)); 1299566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts,th->stage_time)); 1309566063dSJacob Faibussowitsch PetscCall(VecCopy(th->X0,X2)); 1319566063dSJacob Faibussowitsch PetscCall(TSAlpha_SNESSolve(ts,NULL,X2)); 1329566063dSJacob Faibussowitsch PetscCall(TSPostStage(ts,th->stage_time,0,&X2)); 1339566063dSJacob Faibussowitsch PetscCall(TSAdaptCheckStage(ts->adapt,ts,th->stage_time,X2,&stageok)); 134b07a2398SLisandro Dalcin if (!stageok) goto finally; 135b07a2398SLisandro Dalcin 136b07a2398SLisandro Dalcin /* Compute V0 ~ dX/dt at t0 with backward differences */ 1379566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(th->V0)); 1389566063dSJacob Faibussowitsch PetscCall(VecAXPY(th->V0,-3/ts->time_step,X0)); 1399566063dSJacob Faibussowitsch PetscCall(VecAXPY(th->V0,+4/ts->time_step,X1)); 1409566063dSJacob Faibussowitsch PetscCall(VecAXPY(th->V0,-1/ts->time_step,X2)); 141b07a2398SLisandro Dalcin 142b07a2398SLisandro Dalcin /* Rough, lower-order estimate LTE of the initial step */ 1432ffb9264SLisandro Dalcin if (th->vec_lte_work) { 1449566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(th->vec_lte_work)); 1459566063dSJacob Faibussowitsch PetscCall(VecAXPY(th->vec_lte_work,+2,X2)); 1469566063dSJacob Faibussowitsch PetscCall(VecAXPY(th->vec_lte_work,-4,X1)); 1479566063dSJacob Faibussowitsch PetscCall(VecAXPY(th->vec_lte_work,+2,X0)); 148b07a2398SLisandro Dalcin } 149b07a2398SLisandro Dalcin 150b07a2398SLisandro Dalcin finally: 151b07a2398SLisandro Dalcin /* Revert TSAlpha to the initial state (t0,X0) */ 152b07a2398SLisandro Dalcin if (initok) *initok = stageok; 1539566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts,time_step)); 1549566063dSJacob Faibussowitsch PetscCall(TSAlphaSetParams(ts,alpha_m,alpha_f,gamma)); 1559566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol,th->X0)); 156b07a2398SLisandro Dalcin 1579566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X1)); 158b07a2398SLisandro Dalcin PetscFunctionReturn(0); 159b07a2398SLisandro Dalcin } 160b07a2398SLisandro Dalcin 161b07a2398SLisandro Dalcin static PetscErrorCode TSStep_Alpha(TS ts) 162b07a2398SLisandro Dalcin { 163b07a2398SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 164b07a2398SLisandro Dalcin PetscInt rejections = 0; 165b07a2398SLisandro Dalcin PetscBool stageok,accept = PETSC_TRUE; 166b07a2398SLisandro Dalcin PetscReal next_time_step = ts->time_step; 167b07a2398SLisandro Dalcin 168b07a2398SLisandro Dalcin PetscFunctionBegin; 1699566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister(citation,&cited)); 170b07a2398SLisandro Dalcin 171b07a2398SLisandro Dalcin if (!ts->steprollback) { 1729566063dSJacob Faibussowitsch if (th->vec_sol_prev) PetscCall(VecCopy(th->X0,th->vec_sol_prev)); 1739566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol,th->X0)); 1749566063dSJacob Faibussowitsch PetscCall(VecCopy(th->V1,th->V0)); 175b07a2398SLisandro Dalcin } 176b07a2398SLisandro Dalcin 1771566a47fSLisandro Dalcin th->status = TS_STEP_INCOMPLETE; 178b07a2398SLisandro Dalcin while (!ts->reason && th->status != TS_STEP_COMPLETE) { 179b07a2398SLisandro Dalcin 180fecfb714SLisandro Dalcin if (ts->steprestart) { 1819566063dSJacob Faibussowitsch PetscCall(TSAlpha_Restart(ts,&stageok)); 182fecfb714SLisandro Dalcin if (!stageok) goto reject_step; 183b07a2398SLisandro Dalcin } 184b07a2398SLisandro Dalcin 1859566063dSJacob Faibussowitsch PetscCall(TSAlpha_StageTime(ts)); 1869566063dSJacob Faibussowitsch PetscCall(VecCopy(th->X0,th->X1)); 1879566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts,th->stage_time)); 1889566063dSJacob Faibussowitsch PetscCall(TSAlpha_SNESSolve(ts,NULL,th->X1)); 1899566063dSJacob Faibussowitsch PetscCall(TSPostStage(ts,th->stage_time,0,&th->Xa)); 1909566063dSJacob Faibussowitsch PetscCall(TSAdaptCheckStage(ts->adapt,ts,th->stage_time,th->Xa,&stageok)); 191fecfb714SLisandro Dalcin if (!stageok) goto reject_step; 192b07a2398SLisandro Dalcin 1931566a47fSLisandro Dalcin th->status = TS_STEP_PENDING; 1949566063dSJacob Faibussowitsch PetscCall(VecCopy(th->X1,ts->vec_sol)); 1959566063dSJacob Faibussowitsch PetscCall(TSAdaptChoose(ts->adapt,ts,ts->time_step,NULL,&next_time_step,&accept)); 1961566a47fSLisandro Dalcin th->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 197b07a2398SLisandro Dalcin if (!accept) { 1989566063dSJacob Faibussowitsch PetscCall(VecCopy(th->X0,ts->vec_sol)); 199be5899b3SLisandro Dalcin ts->time_step = next_time_step; 200be5899b3SLisandro Dalcin goto reject_step; 201b07a2398SLisandro Dalcin } 202b07a2398SLisandro Dalcin 203b07a2398SLisandro Dalcin ts->ptime += ts->time_step; 204b07a2398SLisandro Dalcin ts->time_step = next_time_step; 205b07a2398SLisandro Dalcin break; 206b07a2398SLisandro Dalcin 207b07a2398SLisandro Dalcin reject_step: 208fecfb714SLisandro Dalcin ts->reject++; accept = PETSC_FALSE; 209b07a2398SLisandro Dalcin if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 210b07a2398SLisandro Dalcin ts->reason = TS_DIVERGED_STEP_REJECTED; 2119566063dSJacob Faibussowitsch PetscCall(PetscInfo(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections)); 212b07a2398SLisandro Dalcin } 2131566a47fSLisandro Dalcin 214b07a2398SLisandro Dalcin } 215b07a2398SLisandro Dalcin PetscFunctionReturn(0); 216b07a2398SLisandro Dalcin } 217b07a2398SLisandro Dalcin 2189808bdc1SLisandro Dalcin static PetscErrorCode TSEvaluateWLTE_Alpha(TS ts,NormType wnormtype,PetscInt *order,PetscReal *wlte) 219b07a2398SLisandro Dalcin { 220b07a2398SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 2219808bdc1SLisandro Dalcin Vec X = th->X1; /* X = solution */ 2221566a47fSLisandro Dalcin Vec Y = th->vec_lte_work; /* Y = X + LTE */ 2237453f775SEmil Constantinescu PetscReal wltea,wlter; 224b07a2398SLisandro Dalcin 225b07a2398SLisandro Dalcin PetscFunctionBegin; 2262ffb9264SLisandro Dalcin if (!th->vec_sol_prev) {*wlte = -1; PetscFunctionReturn(0);} 2272ffb9264SLisandro Dalcin if (!th->vec_lte_work) {*wlte = -1; PetscFunctionReturn(0);} 228fecfb714SLisandro Dalcin if (ts->steprestart) { 229fecfb714SLisandro Dalcin /* th->vec_lte_work is set to the LTE in TSAlpha_Restart() */ 2309566063dSJacob Faibussowitsch PetscCall(VecAXPY(Y,1,X)); 231b07a2398SLisandro Dalcin } else { 232b07a2398SLisandro Dalcin /* Compute LTE using backward differences with non-constant time step */ 233be5899b3SLisandro Dalcin PetscReal h = ts->time_step, h_prev = ts->ptime - ts->ptime_prev; 234be5899b3SLisandro Dalcin PetscReal a = 1 + h_prev/h; 235b07a2398SLisandro Dalcin PetscScalar scal[3]; Vec vecs[3]; 236b07a2398SLisandro Dalcin scal[0] = +1/a; scal[1] = -1/(a-1); scal[2] = +1/(a*(a-1)); 237b07a2398SLisandro Dalcin vecs[0] = th->X1; vecs[1] = th->X0; vecs[2] = th->vec_sol_prev; 2389566063dSJacob Faibussowitsch PetscCall(VecCopy(X,Y)); 2399566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Y,3,scal,vecs)); 240b07a2398SLisandro Dalcin } 2419566063dSJacob Faibussowitsch PetscCall(TSErrorWeightedNorm(ts,X,Y,wnormtype,wlte,&wltea,&wlter)); 2429808bdc1SLisandro Dalcin if (order) *order = 2; 243b07a2398SLisandro Dalcin PetscFunctionReturn(0); 244b07a2398SLisandro Dalcin } 245b07a2398SLisandro Dalcin 246b07a2398SLisandro Dalcin static PetscErrorCode TSRollBack_Alpha(TS ts) 247b07a2398SLisandro Dalcin { 248b07a2398SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 249b07a2398SLisandro Dalcin 250b07a2398SLisandro Dalcin PetscFunctionBegin; 2519566063dSJacob Faibussowitsch PetscCall(VecCopy(th->X0,ts->vec_sol)); 252b07a2398SLisandro Dalcin PetscFunctionReturn(0); 253b07a2398SLisandro Dalcin } 254b07a2398SLisandro Dalcin 255b07a2398SLisandro Dalcin static PetscErrorCode TSInterpolate_Alpha(TS ts,PetscReal t,Vec X) 256b07a2398SLisandro Dalcin { 257b07a2398SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 258b07a2398SLisandro Dalcin PetscReal dt = t - ts->ptime; 259b07a2398SLisandro Dalcin 260b07a2398SLisandro Dalcin PetscFunctionBegin; 2619566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol,X)); 2629566063dSJacob Faibussowitsch PetscCall(VecAXPY(X,th->Gamma*dt,th->V1)); 2639566063dSJacob Faibussowitsch PetscCall(VecAXPY(X,(1-th->Gamma)*dt,th->V0)); 264b07a2398SLisandro Dalcin PetscFunctionReturn(0); 265b07a2398SLisandro Dalcin } 266b07a2398SLisandro Dalcin 267b07a2398SLisandro Dalcin static PetscErrorCode SNESTSFormFunction_Alpha(PETSC_UNUSED SNES snes,Vec X,Vec F,TS ts) 268b07a2398SLisandro Dalcin { 269b07a2398SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 270b07a2398SLisandro Dalcin PetscReal ta = th->stage_time; 271b07a2398SLisandro Dalcin Vec Xa = th->Xa, Va = th->Va; 272b07a2398SLisandro Dalcin 273b07a2398SLisandro Dalcin PetscFunctionBegin; 2749566063dSJacob Faibussowitsch PetscCall(TSAlpha_StageVecs(ts,X)); 275b07a2398SLisandro Dalcin /* F = Function(ta,Xa,Va) */ 2769566063dSJacob Faibussowitsch PetscCall(TSComputeIFunction(ts,ta,Xa,Va,F,PETSC_FALSE)); 2779566063dSJacob Faibussowitsch PetscCall(VecScale(F,th->scale_F)); 278b07a2398SLisandro Dalcin PetscFunctionReturn(0); 279b07a2398SLisandro Dalcin } 280b07a2398SLisandro Dalcin 281b07a2398SLisandro Dalcin static PetscErrorCode SNESTSFormJacobian_Alpha(PETSC_UNUSED SNES snes,PETSC_UNUSED Vec X,Mat J,Mat P,TS ts) 282b07a2398SLisandro Dalcin { 283b07a2398SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 284b07a2398SLisandro Dalcin PetscReal ta = th->stage_time; 285b07a2398SLisandro Dalcin Vec Xa = th->Xa, Va = th->Va; 286b07a2398SLisandro Dalcin PetscReal dVdX = th->shift_V; 287b07a2398SLisandro Dalcin 288b07a2398SLisandro Dalcin PetscFunctionBegin; 289b07a2398SLisandro Dalcin /* J,P = Jacobian(ta,Xa,Va) */ 2909566063dSJacob Faibussowitsch PetscCall(TSComputeIJacobian(ts,ta,Xa,Va,dVdX,J,P,PETSC_FALSE)); 291b07a2398SLisandro Dalcin PetscFunctionReturn(0); 292b07a2398SLisandro Dalcin } 293b07a2398SLisandro Dalcin 294b07a2398SLisandro Dalcin static PetscErrorCode TSReset_Alpha(TS ts) 295b07a2398SLisandro Dalcin { 296b07a2398SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 297b07a2398SLisandro Dalcin 298b07a2398SLisandro Dalcin PetscFunctionBegin; 2999566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->X0)); 3009566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->Xa)); 3019566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->X1)); 3029566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->V0)); 3039566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->Va)); 3049566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->V1)); 3059566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->vec_sol_prev)); 3069566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->vec_lte_work)); 307b07a2398SLisandro Dalcin PetscFunctionReturn(0); 308b07a2398SLisandro Dalcin } 309b07a2398SLisandro Dalcin 310b07a2398SLisandro Dalcin static PetscErrorCode TSDestroy_Alpha(TS ts) 311b07a2398SLisandro Dalcin { 312b07a2398SLisandro Dalcin PetscFunctionBegin; 3139566063dSJacob Faibussowitsch PetscCall(TSReset_Alpha(ts)); 3149566063dSJacob Faibussowitsch PetscCall(PetscFree(ts->data)); 315b07a2398SLisandro Dalcin 3169566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSAlphaSetRadius_C",NULL)); 3179566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSAlphaSetParams_C",NULL)); 3189566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSAlphaGetParams_C",NULL)); 319b07a2398SLisandro Dalcin PetscFunctionReturn(0); 320b07a2398SLisandro Dalcin } 321b07a2398SLisandro Dalcin 322b07a2398SLisandro Dalcin static PetscErrorCode TSSetUp_Alpha(TS ts) 323b07a2398SLisandro Dalcin { 324b07a2398SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 3252ffb9264SLisandro Dalcin PetscBool match; 326b07a2398SLisandro Dalcin 327b07a2398SLisandro Dalcin PetscFunctionBegin; 3289566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol,&th->X0)); 3299566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol,&th->Xa)); 3309566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol,&th->X1)); 3319566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol,&th->V0)); 3329566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol,&th->Va)); 3339566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol,&th->V1)); 3341566a47fSLisandro Dalcin 3359566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts,&ts->adapt)); 3369566063dSJacob Faibussowitsch PetscCall(TSAdaptCandidatesClear(ts->adapt)); 3379566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)ts->adapt,TSADAPTNONE,&match)); 3382ffb9264SLisandro Dalcin if (!match) { 3399566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol,&th->vec_sol_prev)); 3409566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol,&th->vec_lte_work)); 341b07a2398SLisandro Dalcin } 3421566a47fSLisandro Dalcin 3439566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts,&ts->snes)); 344b07a2398SLisandro Dalcin PetscFunctionReturn(0); 345b07a2398SLisandro Dalcin } 346b07a2398SLisandro Dalcin 347b07a2398SLisandro Dalcin static PetscErrorCode TSSetFromOptions_Alpha(PetscOptionItems *PetscOptionsObject,TS ts) 348b07a2398SLisandro Dalcin { 349b07a2398SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 350b07a2398SLisandro Dalcin 351b07a2398SLisandro Dalcin PetscFunctionBegin; 3529566063dSJacob Faibussowitsch PetscCall(PetscOptionsHead(PetscOptionsObject,"Generalized-Alpha ODE solver options")); 353b07a2398SLisandro Dalcin { 354b07a2398SLisandro Dalcin PetscBool flg; 355b07a2398SLisandro Dalcin PetscReal radius = 1; 3569566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_alpha_radius","Spectral radius (high-frequency dissipation)","TSAlphaSetRadius",radius,&radius,&flg)); 3579566063dSJacob Faibussowitsch if (flg) PetscCall(TSAlphaSetRadius(ts,radius)); 3589566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_alpha_alpha_m","Algorithmic parameter alpha_m","TSAlphaSetParams",th->Alpha_m,&th->Alpha_m,NULL)); 3599566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_alpha_alpha_f","Algorithmic parameter alpha_f","TSAlphaSetParams",th->Alpha_f,&th->Alpha_f,NULL)); 3609566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_alpha_gamma","Algorithmic parameter gamma","TSAlphaSetParams",th->Gamma,&th->Gamma,NULL)); 3619566063dSJacob Faibussowitsch PetscCall(TSAlphaSetParams(ts,th->Alpha_m,th->Alpha_f,th->Gamma)); 362b07a2398SLisandro Dalcin } 3639566063dSJacob Faibussowitsch PetscCall(PetscOptionsTail()); 364b07a2398SLisandro Dalcin PetscFunctionReturn(0); 365b07a2398SLisandro Dalcin } 366b07a2398SLisandro Dalcin 367b07a2398SLisandro Dalcin static PetscErrorCode TSView_Alpha(TS ts,PetscViewer viewer) 368b07a2398SLisandro Dalcin { 369b07a2398SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 3709c334d8fSLisandro Dalcin PetscBool iascii; 371b07a2398SLisandro Dalcin 372b07a2398SLisandro Dalcin PetscFunctionBegin; 3739566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 3749c334d8fSLisandro Dalcin if (iascii) { 3759566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Alpha_m=%g, Alpha_f=%g, Gamma=%g\n",(double)th->Alpha_m,(double)th->Alpha_f,(double)th->Gamma)); 3769c334d8fSLisandro Dalcin } 377b07a2398SLisandro Dalcin PetscFunctionReturn(0); 378b07a2398SLisandro Dalcin } 379b07a2398SLisandro Dalcin 380b07a2398SLisandro Dalcin static PetscErrorCode TSAlphaSetRadius_Alpha(TS ts,PetscReal radius) 381b07a2398SLisandro Dalcin { 382b07a2398SLisandro Dalcin PetscReal alpha_m,alpha_f,gamma; 383b07a2398SLisandro Dalcin 384b07a2398SLisandro Dalcin PetscFunctionBegin; 3852c71b3e2SJacob Faibussowitsch PetscCheckFalse(radius < 0 || radius > 1,PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Radius %g not in range [0,1]",(double)radius); 386b07a2398SLisandro Dalcin alpha_m = (PetscReal)0.5*(3-radius)/(1+radius); 387b07a2398SLisandro Dalcin alpha_f = 1/(1+radius); 388b07a2398SLisandro Dalcin gamma = (PetscReal)0.5 + alpha_m - alpha_f; 3899566063dSJacob Faibussowitsch PetscCall(TSAlphaSetParams(ts,alpha_m,alpha_f,gamma)); 390b07a2398SLisandro Dalcin PetscFunctionReturn(0); 391b07a2398SLisandro Dalcin } 392b07a2398SLisandro Dalcin 393b07a2398SLisandro Dalcin static PetscErrorCode TSAlphaSetParams_Alpha(TS ts,PetscReal alpha_m,PetscReal alpha_f,PetscReal gamma) 394b07a2398SLisandro Dalcin { 395b07a2398SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 396b07a2398SLisandro Dalcin PetscReal tol = 100*PETSC_MACHINE_EPSILON; 397b07a2398SLisandro Dalcin PetscReal res = ((PetscReal)0.5 + alpha_m - alpha_f) - gamma; 398b07a2398SLisandro Dalcin 399b07a2398SLisandro Dalcin PetscFunctionBegin; 400b07a2398SLisandro Dalcin th->Alpha_m = alpha_m; 401b07a2398SLisandro Dalcin th->Alpha_f = alpha_f; 402b07a2398SLisandro Dalcin th->Gamma = gamma; 403b07a2398SLisandro Dalcin th->order = (PetscAbsReal(res) < tol) ? 2 : 1; 404b07a2398SLisandro Dalcin PetscFunctionReturn(0); 405b07a2398SLisandro Dalcin } 406b07a2398SLisandro Dalcin 407b07a2398SLisandro Dalcin static PetscErrorCode TSAlphaGetParams_Alpha(TS ts,PetscReal *alpha_m,PetscReal *alpha_f,PetscReal *gamma) 408b07a2398SLisandro Dalcin { 409b07a2398SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 410b07a2398SLisandro Dalcin 411b07a2398SLisandro Dalcin PetscFunctionBegin; 412b07a2398SLisandro Dalcin if (alpha_m) *alpha_m = th->Alpha_m; 413b07a2398SLisandro Dalcin if (alpha_f) *alpha_f = th->Alpha_f; 414b07a2398SLisandro Dalcin if (gamma) *gamma = th->Gamma; 415b07a2398SLisandro Dalcin PetscFunctionReturn(0); 416b07a2398SLisandro Dalcin } 417b07a2398SLisandro Dalcin 418b07a2398SLisandro Dalcin /*MC 419b07a2398SLisandro Dalcin TSALPHA - ODE/DAE solver using the implicit Generalized-Alpha method 420b07a2398SLisandro Dalcin for first-order systems 421b07a2398SLisandro Dalcin 422b07a2398SLisandro Dalcin Level: beginner 423b07a2398SLisandro Dalcin 424b07a2398SLisandro Dalcin References: 425606c0280SSatish Balay + * - K.E. Jansen, C.H. Whiting, G.M. Hulber, "A generalized-alpha 426b07a2398SLisandro Dalcin method for integrating the filtered Navier-Stokes equations with a 427b07a2398SLisandro Dalcin stabilized finite element method", Computer Methods in Applied 428b07a2398SLisandro Dalcin Mechanics and Engineering, 190, 305-319, 2000. 429b07a2398SLisandro Dalcin DOI: 10.1016/S0045-7825(00)00203-6. 430606c0280SSatish Balay - * - J. Chung, G.M.Hubert. "A Time Integration Algorithm for Structural 431b07a2398SLisandro Dalcin Dynamics with Improved Numerical Dissipation: The Generalized-alpha 432b07a2398SLisandro Dalcin Method" ASME Journal of Applied Mechanics, 60, 371:375, 1993. 433b07a2398SLisandro Dalcin 434b07a2398SLisandro Dalcin .seealso: TS, TSCreate(), TSSetType(), TSAlphaSetRadius(), TSAlphaSetParams() 435b07a2398SLisandro Dalcin M*/ 436b07a2398SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSCreate_Alpha(TS ts) 437b07a2398SLisandro Dalcin { 438b07a2398SLisandro Dalcin TS_Alpha *th; 439b07a2398SLisandro Dalcin 440b07a2398SLisandro Dalcin PetscFunctionBegin; 441b07a2398SLisandro Dalcin ts->ops->reset = TSReset_Alpha; 442b07a2398SLisandro Dalcin ts->ops->destroy = TSDestroy_Alpha; 443b07a2398SLisandro Dalcin ts->ops->view = TSView_Alpha; 444b07a2398SLisandro Dalcin ts->ops->setup = TSSetUp_Alpha; 445b07a2398SLisandro Dalcin ts->ops->setfromoptions = TSSetFromOptions_Alpha; 446b07a2398SLisandro Dalcin ts->ops->step = TSStep_Alpha; 4479808bdc1SLisandro Dalcin ts->ops->evaluatewlte = TSEvaluateWLTE_Alpha; 448b07a2398SLisandro Dalcin ts->ops->rollback = TSRollBack_Alpha; 449b07a2398SLisandro Dalcin ts->ops->interpolate = TSInterpolate_Alpha; 450b07a2398SLisandro Dalcin ts->ops->snesfunction = SNESTSFormFunction_Alpha; 451b07a2398SLisandro Dalcin ts->ops->snesjacobian = SNESTSFormJacobian_Alpha; 4522ffb9264SLisandro Dalcin ts->default_adapt_type = TSADAPTNONE; 453b07a2398SLisandro Dalcin 454efd4aadfSBarry Smith ts->usessnes = PETSC_TRUE; 455efd4aadfSBarry Smith 4569566063dSJacob Faibussowitsch PetscCall(PetscNewLog(ts,&th)); 457b07a2398SLisandro Dalcin ts->data = (void*)th; 458b07a2398SLisandro Dalcin 459b07a2398SLisandro Dalcin th->Alpha_m = 0.5; 460b07a2398SLisandro Dalcin th->Alpha_f = 0.5; 461b07a2398SLisandro Dalcin th->Gamma = 0.5; 462b07a2398SLisandro Dalcin th->order = 2; 463b07a2398SLisandro Dalcin 4649566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSAlphaSetRadius_C",TSAlphaSetRadius_Alpha)); 4659566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSAlphaSetParams_C",TSAlphaSetParams_Alpha)); 4669566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSAlphaGetParams_C",TSAlphaGetParams_Alpha)); 467b07a2398SLisandro Dalcin PetscFunctionReturn(0); 468b07a2398SLisandro Dalcin } 469b07a2398SLisandro Dalcin 470b07a2398SLisandro Dalcin /*@ 471b07a2398SLisandro Dalcin TSAlphaSetRadius - sets the desired spectral radius of the method 472b07a2398SLisandro Dalcin (i.e. high-frequency numerical damping) 473b07a2398SLisandro Dalcin 474b07a2398SLisandro Dalcin Logically Collective on TS 475b07a2398SLisandro Dalcin 476b07a2398SLisandro Dalcin The algorithmic parameters \alpha_m and \alpha_f of the 477b07a2398SLisandro Dalcin generalized-\alpha method can be computed in terms of a specified 478b07a2398SLisandro Dalcin spectral radius \rho in [0,1] for infinite time step in order to 479b07a2398SLisandro Dalcin control high-frequency numerical damping: 480b07a2398SLisandro Dalcin \alpha_m = 0.5*(3-\rho)/(1+\rho) 481b07a2398SLisandro Dalcin \alpha_f = 1/(1+\rho) 482b07a2398SLisandro Dalcin 483d8d19677SJose E. Roman Input Parameters: 484b07a2398SLisandro Dalcin + ts - timestepping context 485b07a2398SLisandro Dalcin - radius - the desired spectral radius 486b07a2398SLisandro Dalcin 487b07a2398SLisandro Dalcin Options Database: 48867b8a455SSatish Balay . -ts_alpha_radius <radius> - set alpha radius 489b07a2398SLisandro Dalcin 490b07a2398SLisandro Dalcin Level: intermediate 491b07a2398SLisandro Dalcin 492b07a2398SLisandro Dalcin .seealso: TSAlphaSetParams(), TSAlphaGetParams() 493b07a2398SLisandro Dalcin @*/ 494b07a2398SLisandro Dalcin PetscErrorCode TSAlphaSetRadius(TS ts,PetscReal radius) 495b07a2398SLisandro Dalcin { 496b07a2398SLisandro Dalcin PetscFunctionBegin; 497b07a2398SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 498b07a2398SLisandro Dalcin PetscValidLogicalCollectiveReal(ts,radius,2); 4992c71b3e2SJacob Faibussowitsch PetscCheckFalse(radius < 0 || radius > 1,((PetscObject)ts)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Radius %g not in range [0,1]",(double)radius); 500*cac4c232SBarry Smith PetscTryMethod(ts,"TSAlphaSetRadius_C",(TS,PetscReal),(ts,radius)); 501b07a2398SLisandro Dalcin PetscFunctionReturn(0); 502b07a2398SLisandro Dalcin } 503b07a2398SLisandro Dalcin 504b07a2398SLisandro Dalcin /*@ 505b07a2398SLisandro Dalcin TSAlphaSetParams - sets the algorithmic parameters for TSALPHA 506b07a2398SLisandro Dalcin 507b07a2398SLisandro Dalcin Logically Collective on TS 508b07a2398SLisandro Dalcin 509b07a2398SLisandro Dalcin Second-order accuracy can be obtained so long as: 510b07a2398SLisandro Dalcin \gamma = 0.5 + alpha_m - alpha_f 511b07a2398SLisandro Dalcin 512b07a2398SLisandro Dalcin Unconditional stability requires: 513b07a2398SLisandro Dalcin \alpha_m >= \alpha_f >= 0.5 514b07a2398SLisandro Dalcin 515b07a2398SLisandro Dalcin Backward Euler method is recovered with: 516b07a2398SLisandro Dalcin \alpha_m = \alpha_f = gamma = 1 517b07a2398SLisandro Dalcin 518d8d19677SJose E. Roman Input Parameters: 519b07a2398SLisandro Dalcin + ts - timestepping context 520a5b23f4aSJose E. Roman . \alpha_m - algorithmic parameter 521a5b23f4aSJose E. Roman . \alpha_f - algorithmic parameter 522a5b23f4aSJose E. Roman - \gamma - algorithmic parameter 523b07a2398SLisandro Dalcin 524b07a2398SLisandro Dalcin Options Database: 52567b8a455SSatish Balay + -ts_alpha_alpha_m <alpha_m> - set alpha_m 52667b8a455SSatish Balay . -ts_alpha_alpha_f <alpha_f> - set alpha_f 52767b8a455SSatish Balay - -ts_alpha_gamma <gamma> - set gamma 528b07a2398SLisandro Dalcin 529b07a2398SLisandro Dalcin Note: 530b07a2398SLisandro Dalcin Use of this function is normally only required to hack TSALPHA to 531b07a2398SLisandro Dalcin use a modified integration scheme. Users should call 532b07a2398SLisandro Dalcin TSAlphaSetRadius() to set the desired spectral radius of the methods 533b07a2398SLisandro Dalcin (i.e. high-frequency damping) in order so select optimal values for 534b07a2398SLisandro Dalcin these parameters. 535b07a2398SLisandro Dalcin 536b07a2398SLisandro Dalcin Level: advanced 537b07a2398SLisandro Dalcin 538b07a2398SLisandro Dalcin .seealso: TSAlphaSetRadius(), TSAlphaGetParams() 539b07a2398SLisandro Dalcin @*/ 540b07a2398SLisandro Dalcin PetscErrorCode TSAlphaSetParams(TS ts,PetscReal alpha_m,PetscReal alpha_f,PetscReal gamma) 541b07a2398SLisandro Dalcin { 542b07a2398SLisandro Dalcin PetscFunctionBegin; 543b07a2398SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 544b07a2398SLisandro Dalcin PetscValidLogicalCollectiveReal(ts,alpha_m,2); 545b07a2398SLisandro Dalcin PetscValidLogicalCollectiveReal(ts,alpha_f,3); 546b07a2398SLisandro Dalcin PetscValidLogicalCollectiveReal(ts,gamma,4); 547*cac4c232SBarry Smith PetscTryMethod(ts,"TSAlphaSetParams_C",(TS,PetscReal,PetscReal,PetscReal),(ts,alpha_m,alpha_f,gamma)); 548b07a2398SLisandro Dalcin PetscFunctionReturn(0); 549b07a2398SLisandro Dalcin } 550b07a2398SLisandro Dalcin 551b07a2398SLisandro Dalcin /*@ 552b07a2398SLisandro Dalcin TSAlphaGetParams - gets the algorithmic parameters for TSALPHA 553b07a2398SLisandro Dalcin 554b07a2398SLisandro Dalcin Not Collective 555b07a2398SLisandro Dalcin 556b07a2398SLisandro Dalcin Input Parameter: 557b07a2398SLisandro Dalcin . ts - timestepping context 558b07a2398SLisandro Dalcin 559b07a2398SLisandro Dalcin Output Parameters: 560b07a2398SLisandro Dalcin + \alpha_m - algorithmic parameter 561b07a2398SLisandro Dalcin . \alpha_f - algorithmic parameter 562b07a2398SLisandro Dalcin - \gamma - algorithmic parameter 563b07a2398SLisandro Dalcin 564b07a2398SLisandro Dalcin Note: 565b07a2398SLisandro Dalcin Use of this function is normally only required to hack TSALPHA to 566b07a2398SLisandro Dalcin use a modified integration scheme. Users should call 567b07a2398SLisandro Dalcin TSAlphaSetRadius() to set the high-frequency damping (i.e. spectral 568b07a2398SLisandro Dalcin radius of the method) in order so select optimal values for these 569b07a2398SLisandro Dalcin parameters. 570b07a2398SLisandro Dalcin 571b07a2398SLisandro Dalcin Level: advanced 572b07a2398SLisandro Dalcin 573b07a2398SLisandro Dalcin .seealso: TSAlphaSetRadius(), TSAlphaSetParams() 574b07a2398SLisandro Dalcin @*/ 575b07a2398SLisandro Dalcin PetscErrorCode TSAlphaGetParams(TS ts,PetscReal *alpha_m,PetscReal *alpha_f,PetscReal *gamma) 576b07a2398SLisandro Dalcin { 577b07a2398SLisandro Dalcin PetscFunctionBegin; 578b07a2398SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 579b07a2398SLisandro Dalcin if (alpha_m) PetscValidRealPointer(alpha_m,2); 580b07a2398SLisandro Dalcin if (alpha_f) PetscValidRealPointer(alpha_f,3); 581b07a2398SLisandro Dalcin if (gamma) PetscValidRealPointer(gamma,4); 582*cac4c232SBarry Smith PetscUseMethod(ts,"TSAlphaGetParams_C",(TS,PetscReal*,PetscReal*,PetscReal*),(ts,alpha_m,alpha_f,gamma)); 583b07a2398SLisandro Dalcin PetscFunctionReturn(0); 584b07a2398SLisandro Dalcin } 585