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 PetscErrorCode ierr; 65b07a2398SLisandro Dalcin 66b07a2398SLisandro Dalcin PetscFunctionBegin; 67b07a2398SLisandro Dalcin /* V1 = 1/(Gamma*dT)*(X1-X0) + (1-1/Gamma)*V0 */ 68b07a2398SLisandro Dalcin ierr = VecWAXPY(V1,-1.0,X0,X1);CHKERRQ(ierr); 69b07a2398SLisandro Dalcin ierr = VecAXPBY(V1,1-1/Gamma,1/(Gamma*dt),V0);CHKERRQ(ierr); 70b07a2398SLisandro Dalcin /* Xa = X0 + Alpha_f*(X1-X0) */ 71b07a2398SLisandro Dalcin ierr = VecWAXPY(Xa,-1.0,X0,X1);CHKERRQ(ierr); 72b07a2398SLisandro Dalcin ierr = VecAYPX(Xa,Alpha_f,X0);CHKERRQ(ierr); 73b07a2398SLisandro Dalcin /* Va = V0 + Alpha_m*(V1-V0) */ 74b07a2398SLisandro Dalcin ierr = VecWAXPY(Va,-1.0,V0,V1);CHKERRQ(ierr); 75b07a2398SLisandro Dalcin ierr = VecAYPX(Va,Alpha_m,V0);CHKERRQ(ierr); 76b07a2398SLisandro Dalcin PetscFunctionReturn(0); 77b07a2398SLisandro Dalcin } 78b07a2398SLisandro Dalcin 79b07a2398SLisandro Dalcin static PetscErrorCode TS_SNESSolve(TS ts,Vec b,Vec x) 80b07a2398SLisandro Dalcin { 81b07a2398SLisandro Dalcin PetscInt nits,lits; 82b07a2398SLisandro Dalcin PetscErrorCode ierr; 83b07a2398SLisandro Dalcin 84b07a2398SLisandro Dalcin PetscFunctionBegin; 85b07a2398SLisandro Dalcin ierr = SNESSolve(ts->snes,b,x);CHKERRQ(ierr); 86b07a2398SLisandro Dalcin ierr = SNESGetIterationNumber(ts->snes,&nits);CHKERRQ(ierr); 87b07a2398SLisandro Dalcin ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr); 88b07a2398SLisandro Dalcin ts->snes_its += nits; ts->ksp_its += lits; 89b07a2398SLisandro Dalcin PetscFunctionReturn(0); 90b07a2398SLisandro Dalcin } 91b07a2398SLisandro Dalcin 92b07a2398SLisandro Dalcin /* 93b07a2398SLisandro Dalcin Compute a consistent initial state for the generalized-alpha method. 94b07a2398SLisandro Dalcin - Solve two successive backward Euler steps with halved time step. 95b07a2398SLisandro Dalcin - Compute the initial time derivative using backward differences. 96b07a2398SLisandro Dalcin - If using adaptivity, estimate the LTE of the initial step. 97b07a2398SLisandro Dalcin */ 98fecfb714SLisandro Dalcin static PetscErrorCode TSAlpha_Restart(TS ts,PetscBool *initok) 99b07a2398SLisandro Dalcin { 100b07a2398SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 101b07a2398SLisandro Dalcin PetscReal time_step; 102b07a2398SLisandro Dalcin PetscReal alpha_m,alpha_f,gamma; 103b07a2398SLisandro Dalcin Vec X0 = ts->vec_sol, X1, X2 = th->X1; 104b07a2398SLisandro Dalcin PetscBool stageok; 105b07a2398SLisandro Dalcin PetscErrorCode ierr; 106b07a2398SLisandro Dalcin 107b07a2398SLisandro Dalcin PetscFunctionBegin; 108b07a2398SLisandro Dalcin ierr = VecDuplicate(X0,&X1);CHKERRQ(ierr); 109b07a2398SLisandro Dalcin 110b07a2398SLisandro Dalcin /* Setup backward Euler with halved time step */ 111b07a2398SLisandro Dalcin ierr = TSAlphaGetParams(ts,&alpha_m,&alpha_f,&gamma);CHKERRQ(ierr); 112b07a2398SLisandro Dalcin ierr = TSAlphaSetParams(ts,1,1,1);CHKERRQ(ierr); 113b07a2398SLisandro Dalcin ierr = TSGetTimeStep(ts,&time_step);CHKERRQ(ierr); 114b07a2398SLisandro Dalcin ts->time_step = time_step/2; 115b07a2398SLisandro Dalcin ierr = TSAlpha_StageTime(ts);CHKERRQ(ierr); 116b07a2398SLisandro Dalcin th->stage_time = ts->ptime; 117b07a2398SLisandro Dalcin ierr = VecZeroEntries(th->V0);CHKERRQ(ierr); 118b07a2398SLisandro Dalcin 119b07a2398SLisandro Dalcin /* First BE step, (t0,X0) -> (t1,X1) */ 120b07a2398SLisandro Dalcin th->stage_time += ts->time_step; 121b07a2398SLisandro Dalcin ierr = VecCopy(X0,th->X0);CHKERRQ(ierr); 122b07a2398SLisandro Dalcin ierr = TSPreStage(ts,th->stage_time);CHKERRQ(ierr); 123b07a2398SLisandro Dalcin ierr = VecCopy(th->X0,X1);CHKERRQ(ierr); 124b07a2398SLisandro Dalcin ierr = TS_SNESSolve(ts,NULL,X1);CHKERRQ(ierr); 125b07a2398SLisandro Dalcin ierr = TSPostStage(ts,th->stage_time,0,&X1);CHKERRQ(ierr); 126b07a2398SLisandro Dalcin ierr = TSAdaptCheckStage(ts->adapt,ts,th->stage_time,X1,&stageok);CHKERRQ(ierr); 127b07a2398SLisandro Dalcin if (!stageok) goto finally; 128b07a2398SLisandro Dalcin 129b07a2398SLisandro Dalcin /* Second BE step, (t1,X1) -> (t2,X2) */ 130b07a2398SLisandro Dalcin th->stage_time += ts->time_step; 131b07a2398SLisandro Dalcin ierr = VecCopy(X1,th->X0);CHKERRQ(ierr); 132b07a2398SLisandro Dalcin ierr = TSPreStage(ts,th->stage_time);CHKERRQ(ierr); 133b07a2398SLisandro Dalcin ierr = VecCopy(th->X0,X2);CHKERRQ(ierr); 134b07a2398SLisandro Dalcin ierr = TS_SNESSolve(ts,NULL,X2);CHKERRQ(ierr); 135b07a2398SLisandro Dalcin ierr = TSPostStage(ts,th->stage_time,0,&X2);CHKERRQ(ierr); 136b07a2398SLisandro Dalcin ierr = TSAdaptCheckStage(ts->adapt,ts,th->stage_time,X2,&stageok);CHKERRQ(ierr); 137b07a2398SLisandro Dalcin if (!stageok) goto finally; 138b07a2398SLisandro Dalcin 139b07a2398SLisandro Dalcin /* Compute V0 ~ dX/dt at t0 with backward differences */ 140b07a2398SLisandro Dalcin ierr = VecZeroEntries(th->V0);CHKERRQ(ierr); 141b07a2398SLisandro Dalcin ierr = VecAXPY(th->V0,-3/ts->time_step,X0);CHKERRQ(ierr); 142b07a2398SLisandro Dalcin ierr = VecAXPY(th->V0,+4/ts->time_step,X1);CHKERRQ(ierr); 143b07a2398SLisandro Dalcin ierr = VecAXPY(th->V0,-1/ts->time_step,X2);CHKERRQ(ierr); 144b07a2398SLisandro Dalcin 145b07a2398SLisandro Dalcin /* Rough, lower-order estimate LTE of the initial step */ 1462ffb9264SLisandro Dalcin if (th->vec_lte_work) { 147fecfb714SLisandro Dalcin ierr = VecZeroEntries(th->vec_lte_work);CHKERRQ(ierr); 148fecfb714SLisandro Dalcin ierr = VecAXPY(th->vec_lte_work,+2,X2);CHKERRQ(ierr); 149fecfb714SLisandro Dalcin ierr = VecAXPY(th->vec_lte_work,-4,X1);CHKERRQ(ierr); 150fecfb714SLisandro Dalcin ierr = VecAXPY(th->vec_lte_work,+2,X0);CHKERRQ(ierr); 151b07a2398SLisandro Dalcin } 152b07a2398SLisandro Dalcin 153b07a2398SLisandro Dalcin finally: 154b07a2398SLisandro Dalcin /* Revert TSAlpha to the initial state (t0,X0) */ 155b07a2398SLisandro Dalcin if (initok) *initok = stageok; 156b07a2398SLisandro Dalcin ierr = TSSetTimeStep(ts,time_step);CHKERRQ(ierr); 157b07a2398SLisandro Dalcin ierr = TSAlphaSetParams(ts,alpha_m,alpha_f,gamma);CHKERRQ(ierr); 158b07a2398SLisandro Dalcin ierr = VecCopy(ts->vec_sol,th->X0);CHKERRQ(ierr); 159b07a2398SLisandro Dalcin 160b07a2398SLisandro Dalcin ierr = VecDestroy(&X1);CHKERRQ(ierr); 161b07a2398SLisandro Dalcin PetscFunctionReturn(0); 162b07a2398SLisandro Dalcin } 163b07a2398SLisandro Dalcin 164b07a2398SLisandro Dalcin static PetscErrorCode TSStep_Alpha(TS ts) 165b07a2398SLisandro Dalcin { 166b07a2398SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 167b07a2398SLisandro Dalcin PetscInt rejections = 0; 168b07a2398SLisandro Dalcin PetscBool stageok,accept = PETSC_TRUE; 169b07a2398SLisandro Dalcin PetscReal next_time_step = ts->time_step; 170b07a2398SLisandro Dalcin PetscErrorCode ierr; 171b07a2398SLisandro Dalcin 172b07a2398SLisandro Dalcin PetscFunctionBegin; 173b07a2398SLisandro Dalcin ierr = PetscCitationsRegister(citation,&cited);CHKERRQ(ierr); 174b07a2398SLisandro Dalcin 175b07a2398SLisandro Dalcin if (!ts->steprollback) { 1762ffb9264SLisandro Dalcin if (th->vec_sol_prev) { ierr = VecCopy(th->X0,th->vec_sol_prev);CHKERRQ(ierr); } 177b07a2398SLisandro Dalcin ierr = VecCopy(ts->vec_sol,th->X0);CHKERRQ(ierr); 178b07a2398SLisandro Dalcin ierr = VecCopy(th->V1,th->V0);CHKERRQ(ierr); 179b07a2398SLisandro Dalcin } 180b07a2398SLisandro Dalcin 1811566a47fSLisandro Dalcin th->status = TS_STEP_INCOMPLETE; 182b07a2398SLisandro Dalcin while (!ts->reason && th->status != TS_STEP_COMPLETE) { 183b07a2398SLisandro Dalcin 184fecfb714SLisandro Dalcin if (ts->steprestart) { 185fecfb714SLisandro Dalcin ierr = TSAlpha_Restart(ts,&stageok);CHKERRQ(ierr); 186fecfb714SLisandro Dalcin if (!stageok) goto reject_step; 187b07a2398SLisandro Dalcin } 188b07a2398SLisandro Dalcin 189b07a2398SLisandro Dalcin ierr = TSAlpha_StageTime(ts);CHKERRQ(ierr); 190b07a2398SLisandro Dalcin ierr = VecCopy(th->X0,th->X1);CHKERRQ(ierr); 191b07a2398SLisandro Dalcin ierr = TSPreStage(ts,th->stage_time);CHKERRQ(ierr); 192b07a2398SLisandro Dalcin ierr = TS_SNESSolve(ts,NULL,th->X1);CHKERRQ(ierr); 193be5899b3SLisandro Dalcin ierr = TSPostStage(ts,th->stage_time,0,&th->Xa);CHKERRQ(ierr); 194be5899b3SLisandro Dalcin ierr = TSAdaptCheckStage(ts->adapt,ts,th->stage_time,th->Xa,&stageok);CHKERRQ(ierr); 195fecfb714SLisandro Dalcin if (!stageok) goto reject_step; 196b07a2398SLisandro Dalcin 1971566a47fSLisandro Dalcin th->status = TS_STEP_PENDING; 198b07a2398SLisandro Dalcin ierr = VecCopy(th->X1,ts->vec_sol);CHKERRQ(ierr); 1991566a47fSLisandro Dalcin ierr = TSAdaptChoose(ts->adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr); 2001566a47fSLisandro Dalcin th->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 201b07a2398SLisandro Dalcin if (!accept) { 202b07a2398SLisandro Dalcin ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr); 203be5899b3SLisandro Dalcin ts->time_step = next_time_step; 204be5899b3SLisandro Dalcin goto reject_step; 205b07a2398SLisandro Dalcin } 206b07a2398SLisandro Dalcin 207b07a2398SLisandro Dalcin ts->ptime += ts->time_step; 208b07a2398SLisandro Dalcin ts->time_step = next_time_step; 209b07a2398SLisandro Dalcin break; 210b07a2398SLisandro Dalcin 211b07a2398SLisandro Dalcin reject_step: 212fecfb714SLisandro Dalcin ts->reject++; accept = PETSC_FALSE; 213b07a2398SLisandro Dalcin if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 214b07a2398SLisandro Dalcin ts->reason = TS_DIVERGED_STEP_REJECTED; 215b07a2398SLisandro Dalcin ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr); 216b07a2398SLisandro Dalcin } 2171566a47fSLisandro Dalcin 218b07a2398SLisandro Dalcin } 219b07a2398SLisandro Dalcin PetscFunctionReturn(0); 220b07a2398SLisandro Dalcin } 221b07a2398SLisandro Dalcin 2229808bdc1SLisandro Dalcin static PetscErrorCode TSEvaluateWLTE_Alpha(TS ts,NormType wnormtype,PetscInt *order,PetscReal *wlte) 223b07a2398SLisandro Dalcin { 224b07a2398SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 2259808bdc1SLisandro Dalcin Vec X = th->X1; /* X = solution */ 2261566a47fSLisandro Dalcin Vec Y = th->vec_lte_work; /* Y = X + LTE */ 2277453f775SEmil Constantinescu PetscReal wltea,wlter; 228b07a2398SLisandro Dalcin PetscErrorCode ierr; 229b07a2398SLisandro Dalcin 230b07a2398SLisandro Dalcin PetscFunctionBegin; 2312ffb9264SLisandro Dalcin if (!th->vec_sol_prev) {*wlte = -1; PetscFunctionReturn(0);} 2322ffb9264SLisandro Dalcin if (!th->vec_lte_work) {*wlte = -1; PetscFunctionReturn(0);} 233fecfb714SLisandro Dalcin if (ts->steprestart) { 234fecfb714SLisandro Dalcin /* th->vec_lte_work is set to the LTE in TSAlpha_Restart() */ 235fecfb714SLisandro Dalcin ierr = VecAXPY(Y,1,X);CHKERRQ(ierr); 236b07a2398SLisandro Dalcin } else { 237b07a2398SLisandro Dalcin /* Compute LTE using backward differences with non-constant time step */ 238be5899b3SLisandro Dalcin PetscReal h = ts->time_step, h_prev = ts->ptime - ts->ptime_prev; 239be5899b3SLisandro Dalcin PetscReal a = 1 + h_prev/h; 240b07a2398SLisandro Dalcin PetscScalar scal[3]; Vec vecs[3]; 241b07a2398SLisandro Dalcin scal[0] = +1/a; scal[1] = -1/(a-1); scal[2] = +1/(a*(a-1)); 242b07a2398SLisandro Dalcin vecs[0] = th->X1; vecs[1] = th->X0; vecs[2] = th->vec_sol_prev; 2439808bdc1SLisandro Dalcin ierr = VecCopy(X,Y);CHKERRQ(ierr); 2449808bdc1SLisandro Dalcin ierr = VecMAXPY(Y,3,scal,vecs);CHKERRQ(ierr); 245b07a2398SLisandro Dalcin } 2467453f775SEmil Constantinescu ierr = TSErrorWeightedNorm(ts,X,Y,wnormtype,wlte,&wltea,&wlter);CHKERRQ(ierr); 2479808bdc1SLisandro Dalcin if (order) *order = 2; 248b07a2398SLisandro Dalcin PetscFunctionReturn(0); 249b07a2398SLisandro Dalcin } 250b07a2398SLisandro Dalcin 251b07a2398SLisandro Dalcin static PetscErrorCode TSRollBack_Alpha(TS ts) 252b07a2398SLisandro Dalcin { 253b07a2398SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 254b07a2398SLisandro Dalcin PetscErrorCode ierr; 255b07a2398SLisandro Dalcin 256b07a2398SLisandro Dalcin PetscFunctionBegin; 257b07a2398SLisandro Dalcin ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr); 258b07a2398SLisandro Dalcin PetscFunctionReturn(0); 259b07a2398SLisandro Dalcin } 260b07a2398SLisandro Dalcin 261b07a2398SLisandro Dalcin static PetscErrorCode TSInterpolate_Alpha(TS ts,PetscReal t,Vec X) 262b07a2398SLisandro Dalcin { 263b07a2398SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 264b07a2398SLisandro Dalcin PetscReal dt = t - ts->ptime; 265b07a2398SLisandro Dalcin PetscErrorCode ierr; 266b07a2398SLisandro Dalcin 267b07a2398SLisandro Dalcin PetscFunctionBegin; 268b07a2398SLisandro Dalcin ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 269b07a2398SLisandro Dalcin ierr = VecAXPY(X,th->Gamma*dt,th->V1);CHKERRQ(ierr); 270b07a2398SLisandro Dalcin ierr = VecAXPY(X,(1-th->Gamma)*dt,th->V0);CHKERRQ(ierr); 271b07a2398SLisandro Dalcin PetscFunctionReturn(0); 272b07a2398SLisandro Dalcin } 273b07a2398SLisandro Dalcin 274b07a2398SLisandro Dalcin static PetscErrorCode SNESTSFormFunction_Alpha(PETSC_UNUSED SNES snes,Vec X,Vec F,TS ts) 275b07a2398SLisandro Dalcin { 276b07a2398SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 277b07a2398SLisandro Dalcin PetscReal ta = th->stage_time; 278b07a2398SLisandro Dalcin Vec Xa = th->Xa, Va = th->Va; 279b07a2398SLisandro Dalcin PetscErrorCode ierr; 280b07a2398SLisandro Dalcin 281b07a2398SLisandro Dalcin PetscFunctionBegin; 282b07a2398SLisandro Dalcin ierr = TSAlpha_StageVecs(ts,X);CHKERRQ(ierr); 283b07a2398SLisandro Dalcin /* F = Function(ta,Xa,Va) */ 284b07a2398SLisandro Dalcin ierr = TSComputeIFunction(ts,ta,Xa,Va,F,PETSC_FALSE);CHKERRQ(ierr); 285b07a2398SLisandro Dalcin ierr = VecScale(F,th->scale_F);CHKERRQ(ierr); 286b07a2398SLisandro Dalcin PetscFunctionReturn(0); 287b07a2398SLisandro Dalcin } 288b07a2398SLisandro Dalcin 289b07a2398SLisandro Dalcin static PetscErrorCode SNESTSFormJacobian_Alpha(PETSC_UNUSED SNES snes,PETSC_UNUSED Vec X,Mat J,Mat P,TS ts) 290b07a2398SLisandro Dalcin { 291b07a2398SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 292b07a2398SLisandro Dalcin PetscReal ta = th->stage_time; 293b07a2398SLisandro Dalcin Vec Xa = th->Xa, Va = th->Va; 294b07a2398SLisandro Dalcin PetscReal dVdX = th->shift_V; 295b07a2398SLisandro Dalcin PetscErrorCode ierr; 296b07a2398SLisandro Dalcin 297b07a2398SLisandro Dalcin PetscFunctionBegin; 298b07a2398SLisandro Dalcin /* J,P = Jacobian(ta,Xa,Va) */ 299b07a2398SLisandro Dalcin ierr = TSComputeIJacobian(ts,ta,Xa,Va,dVdX,J,P,PETSC_FALSE);CHKERRQ(ierr); 300b07a2398SLisandro Dalcin PetscFunctionReturn(0); 301b07a2398SLisandro Dalcin } 302b07a2398SLisandro Dalcin 303b07a2398SLisandro Dalcin static PetscErrorCode TSReset_Alpha(TS ts) 304b07a2398SLisandro Dalcin { 305b07a2398SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 306b07a2398SLisandro Dalcin PetscErrorCode ierr; 307b07a2398SLisandro Dalcin 308b07a2398SLisandro Dalcin PetscFunctionBegin; 309b07a2398SLisandro Dalcin ierr = VecDestroy(&th->X0);CHKERRQ(ierr); 310b07a2398SLisandro Dalcin ierr = VecDestroy(&th->Xa);CHKERRQ(ierr); 311b07a2398SLisandro Dalcin ierr = VecDestroy(&th->X1);CHKERRQ(ierr); 312b07a2398SLisandro Dalcin ierr = VecDestroy(&th->V0);CHKERRQ(ierr); 313b07a2398SLisandro Dalcin ierr = VecDestroy(&th->Va);CHKERRQ(ierr); 314b07a2398SLisandro Dalcin ierr = VecDestroy(&th->V1);CHKERRQ(ierr); 315b07a2398SLisandro Dalcin ierr = VecDestroy(&th->vec_sol_prev);CHKERRQ(ierr); 3161566a47fSLisandro Dalcin ierr = VecDestroy(&th->vec_lte_work);CHKERRQ(ierr); 317b07a2398SLisandro Dalcin PetscFunctionReturn(0); 318b07a2398SLisandro Dalcin } 319b07a2398SLisandro Dalcin 320b07a2398SLisandro Dalcin static PetscErrorCode TSDestroy_Alpha(TS ts) 321b07a2398SLisandro Dalcin { 322b07a2398SLisandro Dalcin PetscErrorCode ierr; 323b07a2398SLisandro Dalcin 324b07a2398SLisandro Dalcin PetscFunctionBegin; 325b07a2398SLisandro Dalcin ierr = TSReset_Alpha(ts);CHKERRQ(ierr); 326b07a2398SLisandro Dalcin ierr = PetscFree(ts->data);CHKERRQ(ierr); 327b07a2398SLisandro Dalcin 328b07a2398SLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)ts,"TSAlphaSetRadius_C",NULL);CHKERRQ(ierr); 329b07a2398SLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)ts,"TSAlphaSetParams_C",NULL);CHKERRQ(ierr); 330b07a2398SLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)ts,"TSAlphaGetParams_C",NULL);CHKERRQ(ierr); 331b07a2398SLisandro Dalcin PetscFunctionReturn(0); 332b07a2398SLisandro Dalcin } 333b07a2398SLisandro Dalcin 334b07a2398SLisandro Dalcin static PetscErrorCode TSSetUp_Alpha(TS ts) 335b07a2398SLisandro Dalcin { 336b07a2398SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 3372ffb9264SLisandro Dalcin PetscBool match; 338b07a2398SLisandro Dalcin PetscErrorCode ierr; 339b07a2398SLisandro Dalcin 340b07a2398SLisandro Dalcin PetscFunctionBegin; 341b07a2398SLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->X0);CHKERRQ(ierr); 342b07a2398SLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->Xa);CHKERRQ(ierr); 343b07a2398SLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->X1);CHKERRQ(ierr); 344b07a2398SLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->V0);CHKERRQ(ierr); 345b07a2398SLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->Va);CHKERRQ(ierr); 346b07a2398SLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->V1);CHKERRQ(ierr); 3471566a47fSLisandro Dalcin 348b07a2398SLisandro Dalcin ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr); 3491566a47fSLisandro Dalcin ierr = TSAdaptCandidatesClear(ts->adapt);CHKERRQ(ierr); 3502ffb9264SLisandro Dalcin ierr = PetscObjectTypeCompare((PetscObject)ts->adapt,TSADAPTNONE,&match);CHKERRQ(ierr); 3512ffb9264SLisandro Dalcin if (!match) { 352b07a2398SLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->vec_sol_prev);CHKERRQ(ierr); 3531566a47fSLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->vec_lte_work);CHKERRQ(ierr); 354b07a2398SLisandro Dalcin } 3551566a47fSLisandro Dalcin 356b07a2398SLisandro Dalcin ierr = TSGetSNES(ts,&ts->snes);CHKERRQ(ierr); 357b07a2398SLisandro Dalcin PetscFunctionReturn(0); 358b07a2398SLisandro Dalcin } 359b07a2398SLisandro Dalcin 360b07a2398SLisandro Dalcin static PetscErrorCode TSSetFromOptions_Alpha(PetscOptionItems *PetscOptionsObject,TS ts) 361b07a2398SLisandro Dalcin { 362b07a2398SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 363b07a2398SLisandro Dalcin PetscErrorCode ierr; 364b07a2398SLisandro Dalcin 365b07a2398SLisandro Dalcin PetscFunctionBegin; 366b07a2398SLisandro Dalcin ierr = PetscOptionsHead(PetscOptionsObject,"Generalized-Alpha ODE solver options");CHKERRQ(ierr); 367b07a2398SLisandro Dalcin { 368b07a2398SLisandro Dalcin PetscBool flg; 369b07a2398SLisandro Dalcin PetscReal radius = 1; 370b07a2398SLisandro Dalcin ierr = PetscOptionsReal("-ts_alpha_radius","Spectral radius (high-frequency dissipation)","TSAlphaSetRadius",radius,&radius,&flg);CHKERRQ(ierr); 371b07a2398SLisandro Dalcin if (flg) {ierr = TSAlphaSetRadius(ts,radius);CHKERRQ(ierr);} 372b07a2398SLisandro Dalcin ierr = PetscOptionsReal("-ts_alpha_alpha_m","Algoritmic parameter alpha_m","TSAlphaSetParams",th->Alpha_m,&th->Alpha_m,NULL);CHKERRQ(ierr); 373b07a2398SLisandro Dalcin ierr = PetscOptionsReal("-ts_alpha_alpha_f","Algoritmic parameter alpha_f","TSAlphaSetParams",th->Alpha_f,&th->Alpha_f,NULL);CHKERRQ(ierr); 374b07a2398SLisandro Dalcin ierr = PetscOptionsReal("-ts_alpha_gamma","Algoritmic parameter gamma","TSAlphaSetParams",th->Gamma,&th->Gamma,NULL);CHKERRQ(ierr); 375b07a2398SLisandro Dalcin ierr = TSAlphaSetParams(ts,th->Alpha_m,th->Alpha_f,th->Gamma);CHKERRQ(ierr); 376b07a2398SLisandro Dalcin } 377b07a2398SLisandro Dalcin ierr = PetscOptionsTail();CHKERRQ(ierr); 378b07a2398SLisandro Dalcin PetscFunctionReturn(0); 379b07a2398SLisandro Dalcin } 380b07a2398SLisandro Dalcin 381b07a2398SLisandro Dalcin static PetscErrorCode TSView_Alpha(TS ts,PetscViewer viewer) 382b07a2398SLisandro Dalcin { 383b07a2398SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 3849c334d8fSLisandro Dalcin PetscBool iascii; 385b07a2398SLisandro Dalcin PetscErrorCode ierr; 386b07a2398SLisandro Dalcin 387b07a2398SLisandro Dalcin PetscFunctionBegin; 3889c334d8fSLisandro Dalcin ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 3899c334d8fSLisandro Dalcin if (iascii) { 3909c334d8fSLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," Alpha_m=%g, Alpha_f=%g, Gamma=%g\n",(double)th->Alpha_m,(double)th->Alpha_f,(double)th->Gamma);CHKERRQ(ierr); 3919c334d8fSLisandro Dalcin } 392b07a2398SLisandro Dalcin PetscFunctionReturn(0); 393b07a2398SLisandro Dalcin } 394b07a2398SLisandro Dalcin 395b07a2398SLisandro Dalcin static PetscErrorCode TSAlphaSetRadius_Alpha(TS ts,PetscReal radius) 396b07a2398SLisandro Dalcin { 397b07a2398SLisandro Dalcin PetscReal alpha_m,alpha_f,gamma; 398b07a2398SLisandro Dalcin PetscErrorCode ierr; 399b07a2398SLisandro Dalcin 400b07a2398SLisandro Dalcin PetscFunctionBegin; 401b07a2398SLisandro Dalcin if (radius < 0 || radius > 1) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Radius %g not in range [0,1]",(double)radius); 402b07a2398SLisandro Dalcin alpha_m = (PetscReal)0.5*(3-radius)/(1+radius); 403b07a2398SLisandro Dalcin alpha_f = 1/(1+radius); 404b07a2398SLisandro Dalcin gamma = (PetscReal)0.5 + alpha_m - alpha_f; 405b07a2398SLisandro Dalcin ierr = TSAlphaSetParams(ts,alpha_m,alpha_f,gamma);CHKERRQ(ierr); 406b07a2398SLisandro Dalcin PetscFunctionReturn(0); 407b07a2398SLisandro Dalcin } 408b07a2398SLisandro Dalcin 409b07a2398SLisandro Dalcin static PetscErrorCode TSAlphaSetParams_Alpha(TS ts,PetscReal alpha_m,PetscReal alpha_f,PetscReal gamma) 410b07a2398SLisandro Dalcin { 411b07a2398SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 412b07a2398SLisandro Dalcin PetscReal tol = 100*PETSC_MACHINE_EPSILON; 413b07a2398SLisandro Dalcin PetscReal res = ((PetscReal)0.5 + alpha_m - alpha_f) - gamma; 414b07a2398SLisandro Dalcin 415b07a2398SLisandro Dalcin PetscFunctionBegin; 416b07a2398SLisandro Dalcin th->Alpha_m = alpha_m; 417b07a2398SLisandro Dalcin th->Alpha_f = alpha_f; 418b07a2398SLisandro Dalcin th->Gamma = gamma; 419b07a2398SLisandro Dalcin th->order = (PetscAbsReal(res) < tol) ? 2 : 1; 420b07a2398SLisandro Dalcin PetscFunctionReturn(0); 421b07a2398SLisandro Dalcin } 422b07a2398SLisandro Dalcin 423b07a2398SLisandro Dalcin static PetscErrorCode TSAlphaGetParams_Alpha(TS ts,PetscReal *alpha_m,PetscReal *alpha_f,PetscReal *gamma) 424b07a2398SLisandro Dalcin { 425b07a2398SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 426b07a2398SLisandro Dalcin 427b07a2398SLisandro Dalcin PetscFunctionBegin; 428b07a2398SLisandro Dalcin if (alpha_m) *alpha_m = th->Alpha_m; 429b07a2398SLisandro Dalcin if (alpha_f) *alpha_f = th->Alpha_f; 430b07a2398SLisandro Dalcin if (gamma) *gamma = th->Gamma; 431b07a2398SLisandro Dalcin PetscFunctionReturn(0); 432b07a2398SLisandro Dalcin } 433b07a2398SLisandro Dalcin 434b07a2398SLisandro Dalcin /*MC 435b07a2398SLisandro Dalcin TSALPHA - ODE/DAE solver using the implicit Generalized-Alpha method 436b07a2398SLisandro Dalcin for first-order systems 437b07a2398SLisandro Dalcin 438b07a2398SLisandro Dalcin Level: beginner 439b07a2398SLisandro Dalcin 440b07a2398SLisandro Dalcin References: 441b07a2398SLisandro Dalcin K.E. Jansen, C.H. Whiting, G.M. Hulber, "A generalized-alpha 442b07a2398SLisandro Dalcin method for integrating the filtered Navier-Stokes equations with a 443b07a2398SLisandro Dalcin stabilized finite element method", Computer Methods in Applied 444b07a2398SLisandro Dalcin Mechanics and Engineering, 190, 305-319, 2000. 445b07a2398SLisandro Dalcin DOI: 10.1016/S0045-7825(00)00203-6. 446b07a2398SLisandro Dalcin 447b07a2398SLisandro Dalcin J. Chung, G.M.Hubert. "A Time Integration Algorithm for Structural 448b07a2398SLisandro Dalcin Dynamics with Improved Numerical Dissipation: The Generalized-alpha 449b07a2398SLisandro Dalcin Method" ASME Journal of Applied Mechanics, 60, 371:375, 1993. 450b07a2398SLisandro Dalcin 451b07a2398SLisandro Dalcin .seealso: TS, TSCreate(), TSSetType(), TSAlphaSetRadius(), TSAlphaSetParams() 452b07a2398SLisandro Dalcin M*/ 453b07a2398SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSCreate_Alpha(TS ts) 454b07a2398SLisandro Dalcin { 455b07a2398SLisandro Dalcin TS_Alpha *th; 456b07a2398SLisandro Dalcin PetscErrorCode ierr; 457b07a2398SLisandro Dalcin 458b07a2398SLisandro Dalcin PetscFunctionBegin; 459b07a2398SLisandro Dalcin ts->ops->reset = TSReset_Alpha; 460b07a2398SLisandro Dalcin ts->ops->destroy = TSDestroy_Alpha; 461b07a2398SLisandro Dalcin ts->ops->view = TSView_Alpha; 462b07a2398SLisandro Dalcin ts->ops->setup = TSSetUp_Alpha; 463b07a2398SLisandro Dalcin ts->ops->setfromoptions = TSSetFromOptions_Alpha; 464b07a2398SLisandro Dalcin ts->ops->step = TSStep_Alpha; 4659808bdc1SLisandro Dalcin ts->ops->evaluatewlte = TSEvaluateWLTE_Alpha; 466b07a2398SLisandro Dalcin ts->ops->rollback = TSRollBack_Alpha; 467b07a2398SLisandro Dalcin ts->ops->interpolate = TSInterpolate_Alpha; 468b07a2398SLisandro Dalcin ts->ops->snesfunction = SNESTSFormFunction_Alpha; 469b07a2398SLisandro Dalcin ts->ops->snesjacobian = SNESTSFormJacobian_Alpha; 4702ffb9264SLisandro Dalcin ts->default_adapt_type = TSADAPTNONE; 471b07a2398SLisandro Dalcin 472*efd4aadfSBarry Smith ts->usessnes = PETSC_TRUE; 473*efd4aadfSBarry Smith 474b07a2398SLisandro Dalcin ierr = PetscNewLog(ts,&th);CHKERRQ(ierr); 475b07a2398SLisandro Dalcin ts->data = (void*)th; 476b07a2398SLisandro Dalcin 477b07a2398SLisandro Dalcin th->Alpha_m = 0.5; 478b07a2398SLisandro Dalcin th->Alpha_f = 0.5; 479b07a2398SLisandro Dalcin th->Gamma = 0.5; 480b07a2398SLisandro Dalcin th->order = 2; 481b07a2398SLisandro Dalcin 482b07a2398SLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)ts,"TSAlphaSetRadius_C",TSAlphaSetRadius_Alpha);CHKERRQ(ierr); 483b07a2398SLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)ts,"TSAlphaSetParams_C",TSAlphaSetParams_Alpha);CHKERRQ(ierr); 484b07a2398SLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)ts,"TSAlphaGetParams_C",TSAlphaGetParams_Alpha);CHKERRQ(ierr); 485b07a2398SLisandro Dalcin PetscFunctionReturn(0); 486b07a2398SLisandro Dalcin } 487b07a2398SLisandro Dalcin 488b07a2398SLisandro Dalcin /*@ 489b07a2398SLisandro Dalcin TSAlphaSetRadius - sets the desired spectral radius of the method 490b07a2398SLisandro Dalcin (i.e. high-frequency numerical damping) 491b07a2398SLisandro Dalcin 492b07a2398SLisandro Dalcin Logically Collective on TS 493b07a2398SLisandro Dalcin 494b07a2398SLisandro Dalcin The algorithmic parameters \alpha_m and \alpha_f of the 495b07a2398SLisandro Dalcin generalized-\alpha method can be computed in terms of a specified 496b07a2398SLisandro Dalcin spectral radius \rho in [0,1] for infinite time step in order to 497b07a2398SLisandro Dalcin control high-frequency numerical damping: 498b07a2398SLisandro Dalcin \alpha_m = 0.5*(3-\rho)/(1+\rho) 499b07a2398SLisandro Dalcin \alpha_f = 1/(1+\rho) 500b07a2398SLisandro Dalcin 501b07a2398SLisandro Dalcin Input Parameter: 502b07a2398SLisandro Dalcin + ts - timestepping context 503b07a2398SLisandro Dalcin - radius - the desired spectral radius 504b07a2398SLisandro Dalcin 505b07a2398SLisandro Dalcin Options Database: 506b07a2398SLisandro Dalcin . -ts_alpha_radius <radius> 507b07a2398SLisandro Dalcin 508b07a2398SLisandro Dalcin Level: intermediate 509b07a2398SLisandro Dalcin 510b07a2398SLisandro Dalcin .seealso: TSAlphaSetParams(), TSAlphaGetParams() 511b07a2398SLisandro Dalcin @*/ 512b07a2398SLisandro Dalcin PetscErrorCode TSAlphaSetRadius(TS ts,PetscReal radius) 513b07a2398SLisandro Dalcin { 514b07a2398SLisandro Dalcin PetscErrorCode ierr; 515b07a2398SLisandro Dalcin 516b07a2398SLisandro Dalcin PetscFunctionBegin; 517b07a2398SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 518b07a2398SLisandro Dalcin PetscValidLogicalCollectiveReal(ts,radius,2); 519b07a2398SLisandro Dalcin if (radius < 0 || radius > 1) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Radius %g not in range [0,1]",(double)radius); 520b07a2398SLisandro Dalcin ierr = PetscTryMethod(ts,"TSAlphaSetRadius_C",(TS,PetscReal),(ts,radius));CHKERRQ(ierr); 521b07a2398SLisandro Dalcin PetscFunctionReturn(0); 522b07a2398SLisandro Dalcin } 523b07a2398SLisandro Dalcin 524b07a2398SLisandro Dalcin /*@ 525b07a2398SLisandro Dalcin TSAlphaSetParams - sets the algorithmic parameters for TSALPHA 526b07a2398SLisandro Dalcin 527b07a2398SLisandro Dalcin Logically Collective on TS 528b07a2398SLisandro Dalcin 529b07a2398SLisandro Dalcin Second-order accuracy can be obtained so long as: 530b07a2398SLisandro Dalcin \gamma = 0.5 + alpha_m - alpha_f 531b07a2398SLisandro Dalcin 532b07a2398SLisandro Dalcin Unconditional stability requires: 533b07a2398SLisandro Dalcin \alpha_m >= \alpha_f >= 0.5 534b07a2398SLisandro Dalcin 535b07a2398SLisandro Dalcin Backward Euler method is recovered with: 536b07a2398SLisandro Dalcin \alpha_m = \alpha_f = gamma = 1 537b07a2398SLisandro Dalcin 538b07a2398SLisandro Dalcin Input Parameter: 539b07a2398SLisandro Dalcin + ts - timestepping context 540b07a2398SLisandro Dalcin . \alpha_m - algorithmic paramenter 541b07a2398SLisandro Dalcin . \alpha_f - algorithmic paramenter 542b07a2398SLisandro Dalcin - \gamma - algorithmic paramenter 543b07a2398SLisandro Dalcin 544b07a2398SLisandro Dalcin Options Database: 545b07a2398SLisandro Dalcin + -ts_alpha_alpha_m <alpha_m> 546b07a2398SLisandro Dalcin . -ts_alpha_alpha_f <alpha_f> 547b07a2398SLisandro Dalcin - -ts_alpha_gamma <gamma> 548b07a2398SLisandro Dalcin 549b07a2398SLisandro Dalcin Note: 550b07a2398SLisandro Dalcin Use of this function is normally only required to hack TSALPHA to 551b07a2398SLisandro Dalcin use a modified integration scheme. Users should call 552b07a2398SLisandro Dalcin TSAlphaSetRadius() to set the desired spectral radius of the methods 553b07a2398SLisandro Dalcin (i.e. high-frequency damping) in order so select optimal values for 554b07a2398SLisandro Dalcin these parameters. 555b07a2398SLisandro Dalcin 556b07a2398SLisandro Dalcin Level: advanced 557b07a2398SLisandro Dalcin 558b07a2398SLisandro Dalcin .seealso: TSAlphaSetRadius(), TSAlphaGetParams() 559b07a2398SLisandro Dalcin @*/ 560b07a2398SLisandro Dalcin PetscErrorCode TSAlphaSetParams(TS ts,PetscReal alpha_m,PetscReal alpha_f,PetscReal gamma) 561b07a2398SLisandro Dalcin { 562b07a2398SLisandro Dalcin PetscErrorCode ierr; 563b07a2398SLisandro Dalcin 564b07a2398SLisandro Dalcin PetscFunctionBegin; 565b07a2398SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 566b07a2398SLisandro Dalcin PetscValidLogicalCollectiveReal(ts,alpha_m,2); 567b07a2398SLisandro Dalcin PetscValidLogicalCollectiveReal(ts,alpha_f,3); 568b07a2398SLisandro Dalcin PetscValidLogicalCollectiveReal(ts,gamma,4); 569b07a2398SLisandro Dalcin ierr = PetscTryMethod(ts,"TSAlphaSetParams_C",(TS,PetscReal,PetscReal,PetscReal),(ts,alpha_m,alpha_f,gamma));CHKERRQ(ierr); 570b07a2398SLisandro Dalcin PetscFunctionReturn(0); 571b07a2398SLisandro Dalcin } 572b07a2398SLisandro Dalcin 573b07a2398SLisandro Dalcin /*@ 574b07a2398SLisandro Dalcin TSAlphaGetParams - gets the algorithmic parameters for TSALPHA 575b07a2398SLisandro Dalcin 576b07a2398SLisandro Dalcin Not Collective 577b07a2398SLisandro Dalcin 578b07a2398SLisandro Dalcin Input Parameter: 579b07a2398SLisandro Dalcin . ts - timestepping context 580b07a2398SLisandro Dalcin 581b07a2398SLisandro Dalcin Output Parameters: 582b07a2398SLisandro Dalcin + \alpha_m - algorithmic parameter 583b07a2398SLisandro Dalcin . \alpha_f - algorithmic parameter 584b07a2398SLisandro Dalcin - \gamma - algorithmic parameter 585b07a2398SLisandro Dalcin 586b07a2398SLisandro Dalcin Note: 587b07a2398SLisandro Dalcin Use of this function is normally only required to hack TSALPHA to 588b07a2398SLisandro Dalcin use a modified integration scheme. Users should call 589b07a2398SLisandro Dalcin TSAlphaSetRadius() to set the high-frequency damping (i.e. spectral 590b07a2398SLisandro Dalcin radius of the method) in order so select optimal values for these 591b07a2398SLisandro Dalcin parameters. 592b07a2398SLisandro Dalcin 593b07a2398SLisandro Dalcin Level: advanced 594b07a2398SLisandro Dalcin 595b07a2398SLisandro Dalcin .seealso: TSAlphaSetRadius(), TSAlphaSetParams() 596b07a2398SLisandro Dalcin @*/ 597b07a2398SLisandro Dalcin PetscErrorCode TSAlphaGetParams(TS ts,PetscReal *alpha_m,PetscReal *alpha_f,PetscReal *gamma) 598b07a2398SLisandro Dalcin { 599b07a2398SLisandro Dalcin PetscErrorCode ierr; 600b07a2398SLisandro Dalcin 601b07a2398SLisandro Dalcin PetscFunctionBegin; 602b07a2398SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 603b07a2398SLisandro Dalcin if (alpha_m) PetscValidRealPointer(alpha_m,2); 604b07a2398SLisandro Dalcin if (alpha_f) PetscValidRealPointer(alpha_f,3); 605b07a2398SLisandro Dalcin if (gamma) PetscValidRealPointer(gamma,4); 606b07a2398SLisandro Dalcin ierr = PetscUseMethod(ts,"TSAlphaGetParams_C",(TS,PetscReal*,PetscReal*,PetscReal*),(ts,alpha_m,alpha_f,gamma));CHKERRQ(ierr); 607b07a2398SLisandro Dalcin PetscFunctionReturn(0); 608b07a2398SLisandro Dalcin } 609