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; 8818efac9SLisandro Dalcin static const char citation[] = 9818efac9SLisandro Dalcin "@article{Chung1993,\n" 10818efac9SLisandro Dalcin " title = {A Time Integration Algorithm for Structural Dynamics with Improved Numerical Dissipation: The Generalized-$\\alpha$ Method},\n" 11818efac9SLisandro Dalcin " author = {J. Chung, G. M. Hubert},\n" 12818efac9SLisandro Dalcin " journal = {ASME Journal of Applied Mechanics},\n" 13818efac9SLisandro Dalcin " volume = {60},\n" 14818efac9SLisandro Dalcin " number = {2},\n" 15818efac9SLisandro Dalcin " pages = {371--375},\n" 16818efac9SLisandro Dalcin " year = {1993},\n" 17818efac9SLisandro Dalcin " issn = {0021-8936},\n" 18818efac9SLisandro Dalcin " doi = {http://dx.doi.org/10.1115/1.2900803}\n}\n"; 19818efac9SLisandro Dalcin 20818efac9SLisandro Dalcin typedef struct { 21818efac9SLisandro Dalcin PetscReal stage_time; 22818efac9SLisandro Dalcin PetscReal shift_V; 23818efac9SLisandro Dalcin PetscReal shift_A; 24818efac9SLisandro Dalcin PetscReal scale_F; 25818efac9SLisandro Dalcin Vec X0,Xa,X1; 26818efac9SLisandro Dalcin Vec V0,Va,V1; 27818efac9SLisandro Dalcin Vec A0,Aa,A1; 28818efac9SLisandro Dalcin 29818efac9SLisandro Dalcin Vec vec_dot; 30818efac9SLisandro Dalcin 31818efac9SLisandro Dalcin PetscReal Alpha_m; 32818efac9SLisandro Dalcin PetscReal Alpha_f; 33818efac9SLisandro Dalcin PetscReal Gamma; 34818efac9SLisandro Dalcin PetscReal Beta; 35818efac9SLisandro Dalcin PetscInt order; 36818efac9SLisandro Dalcin 37818efac9SLisandro Dalcin Vec vec_sol_prev; 38818efac9SLisandro Dalcin Vec vec_dot_prev; 39818efac9SLisandro Dalcin Vec vec_lte_work[2]; 40818efac9SLisandro Dalcin 41818efac9SLisandro Dalcin TSStepStatus status; 42818efac9SLisandro Dalcin } TS_Alpha; 43818efac9SLisandro Dalcin 44818efac9SLisandro Dalcin static PetscErrorCode TSAlpha_StageTime(TS ts) 45818efac9SLisandro Dalcin { 46818efac9SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 47818efac9SLisandro Dalcin PetscReal t = ts->ptime; 48818efac9SLisandro Dalcin PetscReal dt = ts->time_step; 49818efac9SLisandro Dalcin PetscReal Alpha_m = th->Alpha_m; 50818efac9SLisandro Dalcin PetscReal Alpha_f = th->Alpha_f; 51818efac9SLisandro Dalcin PetscReal Gamma = th->Gamma; 52818efac9SLisandro Dalcin PetscReal Beta = th->Beta; 53818efac9SLisandro Dalcin 54818efac9SLisandro Dalcin PetscFunctionBegin; 55818efac9SLisandro Dalcin th->stage_time = t + Alpha_f*dt; 56818efac9SLisandro Dalcin th->shift_V = Gamma/(dt*Beta); 57818efac9SLisandro Dalcin th->shift_A = Alpha_m/(Alpha_f*dt*dt*Beta); 58818efac9SLisandro Dalcin th->scale_F = 1/Alpha_f; 59818efac9SLisandro Dalcin PetscFunctionReturn(0); 60818efac9SLisandro Dalcin } 61818efac9SLisandro Dalcin 62818efac9SLisandro Dalcin static PetscErrorCode TSAlpha_StageVecs(TS ts,Vec X) 63818efac9SLisandro Dalcin { 64818efac9SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 65818efac9SLisandro Dalcin Vec X1 = X, V1 = th->V1, A1 = th->A1; 66818efac9SLisandro Dalcin Vec Xa = th->Xa, Va = th->Va, Aa = th->Aa; 67818efac9SLisandro Dalcin Vec X0 = th->X0, V0 = th->V0, A0 = th->A0; 68818efac9SLisandro Dalcin PetscReal dt = ts->time_step; 69818efac9SLisandro Dalcin PetscReal Alpha_m = th->Alpha_m; 70818efac9SLisandro Dalcin PetscReal Alpha_f = th->Alpha_f; 71818efac9SLisandro Dalcin PetscReal Gamma = th->Gamma; 72818efac9SLisandro Dalcin PetscReal Beta = th->Beta; 73818efac9SLisandro Dalcin 74818efac9SLisandro Dalcin PetscFunctionBegin; 75818efac9SLisandro Dalcin /* A1 = ... */ 769566063dSJacob Faibussowitsch PetscCall(VecWAXPY(A1,-1.0,X0,X1)); 779566063dSJacob Faibussowitsch PetscCall(VecAXPY (A1,-dt,V0)); 789566063dSJacob Faibussowitsch PetscCall(VecAXPBY(A1,-(1-2*Beta)/(2*Beta),1/(dt*dt*Beta),A0)); 79818efac9SLisandro Dalcin /* V1 = ... */ 809566063dSJacob Faibussowitsch PetscCall(VecWAXPY(V1,(1.0-Gamma)/Gamma,A0,A1)); 819566063dSJacob Faibussowitsch PetscCall(VecAYPX (V1,dt*Gamma,V0)); 82818efac9SLisandro Dalcin /* Xa = X0 + Alpha_f*(X1-X0) */ 839566063dSJacob Faibussowitsch PetscCall(VecWAXPY(Xa,-1.0,X0,X1)); 849566063dSJacob Faibussowitsch PetscCall(VecAYPX (Xa,Alpha_f,X0)); 85818efac9SLisandro Dalcin /* Va = V0 + Alpha_f*(V1-V0) */ 869566063dSJacob Faibussowitsch PetscCall(VecWAXPY(Va,-1.0,V0,V1)); 879566063dSJacob Faibussowitsch PetscCall(VecAYPX (Va,Alpha_f,V0)); 88818efac9SLisandro Dalcin /* Aa = A0 + Alpha_m*(A1-A0) */ 899566063dSJacob Faibussowitsch PetscCall(VecWAXPY(Aa,-1.0,A0,A1)); 909566063dSJacob Faibussowitsch PetscCall(VecAYPX (Aa,Alpha_m,A0)); 91818efac9SLisandro Dalcin PetscFunctionReturn(0); 92818efac9SLisandro Dalcin } 93818efac9SLisandro Dalcin 94470880abSPatrick Sanan static PetscErrorCode TSAlpha_SNESSolve(TS ts,Vec b,Vec x) 95818efac9SLisandro Dalcin { 96818efac9SLisandro Dalcin PetscInt nits,lits; 97818efac9SLisandro Dalcin 98818efac9SLisandro Dalcin PetscFunctionBegin; 999566063dSJacob Faibussowitsch PetscCall(SNESSolve(ts->snes,b,x)); 1009566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(ts->snes,&nits)); 1019566063dSJacob Faibussowitsch PetscCall(SNESGetLinearSolveIterations(ts->snes,&lits)); 102818efac9SLisandro Dalcin ts->snes_its += nits; ts->ksp_its += lits; 103818efac9SLisandro Dalcin PetscFunctionReturn(0); 104818efac9SLisandro Dalcin } 105818efac9SLisandro Dalcin 106818efac9SLisandro Dalcin /* 107818efac9SLisandro Dalcin Compute a consistent initial state for the generalized-alpha method. 108818efac9SLisandro Dalcin - Solve two successive backward Euler steps with halved time step. 109818efac9SLisandro Dalcin - Compute the initial second time derivative using backward differences. 110818efac9SLisandro Dalcin - If using adaptivity, estimate the LTE of the initial step. 111818efac9SLisandro Dalcin */ 112818efac9SLisandro Dalcin static PetscErrorCode TSAlpha_Restart(TS ts,PetscBool *initok) 113818efac9SLisandro Dalcin { 114818efac9SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 115818efac9SLisandro Dalcin PetscReal time_step; 116818efac9SLisandro Dalcin PetscReal alpha_m,alpha_f,gamma,beta; 117818efac9SLisandro Dalcin Vec X0 = ts->vec_sol, X1, X2 = th->X1; 118818efac9SLisandro Dalcin Vec V0 = ts->vec_dot, V1, V2 = th->V1; 119818efac9SLisandro Dalcin PetscBool stageok; 120818efac9SLisandro Dalcin 121818efac9SLisandro Dalcin PetscFunctionBegin; 1229566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X0,&X1)); 1239566063dSJacob Faibussowitsch PetscCall(VecDuplicate(V0,&V1)); 124818efac9SLisandro Dalcin 125818efac9SLisandro Dalcin /* Setup backward Euler with halved time step */ 1269566063dSJacob Faibussowitsch PetscCall(TSAlpha2GetParams(ts,&alpha_m,&alpha_f,&gamma,&beta)); 1279566063dSJacob Faibussowitsch PetscCall(TSAlpha2SetParams(ts,1,1,1,0.5)); 1289566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts,&time_step)); 129818efac9SLisandro Dalcin ts->time_step = time_step/2; 1309566063dSJacob Faibussowitsch PetscCall(TSAlpha_StageTime(ts)); 131818efac9SLisandro Dalcin th->stage_time = ts->ptime; 1329566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(th->A0)); 133818efac9SLisandro Dalcin 134818efac9SLisandro Dalcin /* First BE step, (t0,X0,V0) -> (t1,X1,V1) */ 135818efac9SLisandro Dalcin th->stage_time += ts->time_step; 1369566063dSJacob Faibussowitsch PetscCall(VecCopy(X0,th->X0)); 1379566063dSJacob Faibussowitsch PetscCall(VecCopy(V0,th->V0)); 1389566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts,th->stage_time)); 1399566063dSJacob Faibussowitsch PetscCall(VecCopy(th->X0,X1)); 1409566063dSJacob Faibussowitsch PetscCall(TSAlpha_SNESSolve(ts,NULL,X1)); 1419566063dSJacob Faibussowitsch PetscCall(VecCopy(th->V1,V1)); 1429566063dSJacob Faibussowitsch PetscCall(TSPostStage(ts,th->stage_time,0,&X1)); 1439566063dSJacob Faibussowitsch PetscCall(TSAdaptCheckStage(ts->adapt,ts,th->stage_time,X1,&stageok)); 144818efac9SLisandro Dalcin if (!stageok) goto finally; 145818efac9SLisandro Dalcin 146818efac9SLisandro Dalcin /* Second BE step, (t1,X1,V1) -> (t2,X2,V2) */ 147818efac9SLisandro Dalcin th->stage_time += ts->time_step; 1489566063dSJacob Faibussowitsch PetscCall(VecCopy(X1,th->X0)); 1499566063dSJacob Faibussowitsch PetscCall(VecCopy(V1,th->V0)); 1509566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts,th->stage_time)); 1519566063dSJacob Faibussowitsch PetscCall(VecCopy(th->X0,X2)); 1529566063dSJacob Faibussowitsch PetscCall(TSAlpha_SNESSolve(ts,NULL,X2)); 1539566063dSJacob Faibussowitsch PetscCall(VecCopy(th->V1,V2)); 1549566063dSJacob Faibussowitsch PetscCall(TSPostStage(ts,th->stage_time,0,&X2)); 1559566063dSJacob Faibussowitsch PetscCall(TSAdaptCheckStage(ts->adapt,ts,th->stage_time,X1,&stageok)); 156818efac9SLisandro Dalcin if (!stageok) goto finally; 157818efac9SLisandro Dalcin 158818efac9SLisandro Dalcin /* Compute A0 ~ dV/dt at t0 with backward differences */ 1599566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(th->A0)); 1609566063dSJacob Faibussowitsch PetscCall(VecAXPY(th->A0,-3/ts->time_step,V0)); 1619566063dSJacob Faibussowitsch PetscCall(VecAXPY(th->A0,+4/ts->time_step,V1)); 1629566063dSJacob Faibussowitsch PetscCall(VecAXPY(th->A0,-1/ts->time_step,V2)); 163818efac9SLisandro Dalcin 164818efac9SLisandro Dalcin /* Rough, lower-order estimate LTE of the initial step */ 1652ffb9264SLisandro Dalcin if (th->vec_lte_work[0]) { 1669566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(th->vec_lte_work[0])); 1679566063dSJacob Faibussowitsch PetscCall(VecAXPY(th->vec_lte_work[0],+2,X2)); 1689566063dSJacob Faibussowitsch PetscCall(VecAXPY(th->vec_lte_work[0],-4,X1)); 1699566063dSJacob Faibussowitsch PetscCall(VecAXPY(th->vec_lte_work[0],+2,X0)); 170818efac9SLisandro Dalcin } 1712ffb9264SLisandro Dalcin if (th->vec_lte_work[1]) { 1729566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(th->vec_lte_work[1])); 1739566063dSJacob Faibussowitsch PetscCall(VecAXPY(th->vec_lte_work[1],+2,V2)); 1749566063dSJacob Faibussowitsch PetscCall(VecAXPY(th->vec_lte_work[1],-4,V1)); 1759566063dSJacob Faibussowitsch PetscCall(VecAXPY(th->vec_lte_work[1],+2,V0)); 176818efac9SLisandro Dalcin } 177818efac9SLisandro Dalcin 178818efac9SLisandro Dalcin finally: 179818efac9SLisandro Dalcin /* Revert TSAlpha to the initial state (t0,X0,V0) */ 180818efac9SLisandro Dalcin if (initok) *initok = stageok; 1819566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts,time_step)); 1829566063dSJacob Faibussowitsch PetscCall(TSAlpha2SetParams(ts,alpha_m,alpha_f,gamma,beta)); 1839566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol,th->X0)); 1849566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_dot,th->V0)); 185818efac9SLisandro Dalcin 1869566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X1)); 1879566063dSJacob Faibussowitsch PetscCall(VecDestroy(&V1)); 188818efac9SLisandro Dalcin PetscFunctionReturn(0); 189818efac9SLisandro Dalcin } 190818efac9SLisandro Dalcin 191818efac9SLisandro Dalcin static PetscErrorCode TSStep_Alpha(TS ts) 192818efac9SLisandro Dalcin { 193818efac9SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 194818efac9SLisandro Dalcin PetscInt rejections = 0; 195818efac9SLisandro Dalcin PetscBool stageok,accept = PETSC_TRUE; 196818efac9SLisandro Dalcin PetscReal next_time_step = ts->time_step; 197818efac9SLisandro Dalcin 198818efac9SLisandro Dalcin PetscFunctionBegin; 1999566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister(citation,&cited)); 200818efac9SLisandro Dalcin 201818efac9SLisandro Dalcin if (!ts->steprollback) { 2029566063dSJacob Faibussowitsch if (th->vec_sol_prev) PetscCall(VecCopy(th->X0,th->vec_sol_prev)); 2039566063dSJacob Faibussowitsch if (th->vec_dot_prev) PetscCall(VecCopy(th->V0,th->vec_dot_prev)); 2049566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol,th->X0)); 2059566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_dot,th->V0)); 2069566063dSJacob Faibussowitsch PetscCall(VecCopy(th->A1,th->A0)); 207818efac9SLisandro Dalcin } 208818efac9SLisandro Dalcin 209818efac9SLisandro Dalcin th->status = TS_STEP_INCOMPLETE; 210818efac9SLisandro Dalcin while (!ts->reason && th->status != TS_STEP_COMPLETE) { 211818efac9SLisandro Dalcin 212818efac9SLisandro Dalcin if (ts->steprestart) { 2139566063dSJacob Faibussowitsch PetscCall(TSAlpha_Restart(ts,&stageok)); 214818efac9SLisandro Dalcin if (!stageok) goto reject_step; 215818efac9SLisandro Dalcin } 216818efac9SLisandro Dalcin 2179566063dSJacob Faibussowitsch PetscCall(TSAlpha_StageTime(ts)); 2189566063dSJacob Faibussowitsch PetscCall(VecCopy(th->X0,th->X1)); 2199566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts,th->stage_time)); 2209566063dSJacob Faibussowitsch PetscCall(TSAlpha_SNESSolve(ts,NULL,th->X1)); 2219566063dSJacob Faibussowitsch PetscCall(TSPostStage(ts,th->stage_time,0,&th->Xa)); 2229566063dSJacob Faibussowitsch PetscCall(TSAdaptCheckStage(ts->adapt,ts,th->stage_time,th->Xa,&stageok)); 223818efac9SLisandro Dalcin if (!stageok) goto reject_step; 224818efac9SLisandro Dalcin 225818efac9SLisandro Dalcin th->status = TS_STEP_PENDING; 2269566063dSJacob Faibussowitsch PetscCall(VecCopy(th->X1,ts->vec_sol)); 2279566063dSJacob Faibussowitsch PetscCall(VecCopy(th->V1,ts->vec_dot)); 2289566063dSJacob Faibussowitsch PetscCall(TSAdaptChoose(ts->adapt,ts,ts->time_step,NULL,&next_time_step,&accept)); 229818efac9SLisandro Dalcin th->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 230818efac9SLisandro Dalcin if (!accept) { 2319566063dSJacob Faibussowitsch PetscCall(VecCopy(th->X0,ts->vec_sol)); 2329566063dSJacob Faibussowitsch PetscCall(VecCopy(th->V0,ts->vec_dot)); 233818efac9SLisandro Dalcin ts->time_step = next_time_step; 234818efac9SLisandro Dalcin goto reject_step; 235818efac9SLisandro Dalcin } 236818efac9SLisandro Dalcin 237818efac9SLisandro Dalcin ts->ptime += ts->time_step; 238818efac9SLisandro Dalcin ts->time_step = next_time_step; 239818efac9SLisandro Dalcin break; 240818efac9SLisandro Dalcin 241818efac9SLisandro Dalcin reject_step: 242818efac9SLisandro Dalcin ts->reject++; accept = PETSC_FALSE; 243818efac9SLisandro Dalcin if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 244818efac9SLisandro Dalcin ts->reason = TS_DIVERGED_STEP_REJECTED; 2459566063dSJacob Faibussowitsch PetscCall(PetscInfo(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections)); 246818efac9SLisandro Dalcin } 247818efac9SLisandro Dalcin 248818efac9SLisandro Dalcin } 249818efac9SLisandro Dalcin PetscFunctionReturn(0); 250818efac9SLisandro Dalcin } 251818efac9SLisandro Dalcin 252818efac9SLisandro Dalcin static PetscErrorCode TSEvaluateWLTE_Alpha(TS ts,NormType wnormtype,PetscInt *order,PetscReal *wlte) 253818efac9SLisandro Dalcin { 254818efac9SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 255818efac9SLisandro Dalcin Vec X = th->X1; /* X = solution */ 256818efac9SLisandro Dalcin Vec V = th->V1; /* V = solution */ 257818efac9SLisandro Dalcin Vec Y = th->vec_lte_work[0]; /* Y = X + LTE */ 258818efac9SLisandro Dalcin Vec Z = th->vec_lte_work[1]; /* Z = V + LTE */ 2597453f775SEmil Constantinescu PetscReal enormX,enormV,enormXa,enormVa,enormXr,enormVr; 260818efac9SLisandro Dalcin 261818efac9SLisandro Dalcin PetscFunctionBegin; 2622ffb9264SLisandro Dalcin if (!th->vec_sol_prev) {*wlte = -1; PetscFunctionReturn(0);} 2632ffb9264SLisandro Dalcin if (!th->vec_dot_prev) {*wlte = -1; PetscFunctionReturn(0);} 2642ffb9264SLisandro Dalcin if (!th->vec_lte_work[0]) {*wlte = -1; PetscFunctionReturn(0);} 2652ffb9264SLisandro Dalcin if (!th->vec_lte_work[1]) {*wlte = -1; PetscFunctionReturn(0);} 266818efac9SLisandro Dalcin if (ts->steprestart) { 2672ffb9264SLisandro Dalcin /* th->vec_lte_prev is set to the LTE in TSAlpha_Restart() */ 2689566063dSJacob Faibussowitsch PetscCall(VecAXPY(Y,1,X)); 2699566063dSJacob Faibussowitsch PetscCall(VecAXPY(Z,1,V)); 270818efac9SLisandro Dalcin } else { 271818efac9SLisandro Dalcin /* Compute LTE using backward differences with non-constant time step */ 272818efac9SLisandro Dalcin PetscReal h = ts->time_step, h_prev = ts->ptime - ts->ptime_prev; 273818efac9SLisandro Dalcin PetscReal a = 1 + h_prev/h; 274818efac9SLisandro Dalcin PetscScalar scal[3]; Vec vecX[3],vecV[3]; 275818efac9SLisandro Dalcin scal[0] = +1/a; scal[1] = -1/(a-1); scal[2] = +1/(a*(a-1)); 276818efac9SLisandro Dalcin vecX[0] = th->X1; vecX[1] = th->X0; vecX[2] = th->vec_sol_prev; 277818efac9SLisandro Dalcin vecV[0] = th->V1; vecV[1] = th->V0; vecV[2] = th->vec_dot_prev; 2789566063dSJacob Faibussowitsch PetscCall(VecCopy(X,Y)); 2799566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Y,3,scal,vecX)); 2809566063dSJacob Faibussowitsch PetscCall(VecCopy(V,Z)); 2819566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Z,3,scal,vecV)); 282818efac9SLisandro Dalcin } 283818efac9SLisandro Dalcin /* XXX ts->atol and ts->vatol are not appropriate for computing enormV */ 2849566063dSJacob Faibussowitsch PetscCall(TSErrorWeightedNorm(ts,X,Y,wnormtype,&enormX,&enormXa,&enormXr)); 2859566063dSJacob Faibussowitsch PetscCall(TSErrorWeightedNorm(ts,V,Z,wnormtype,&enormV,&enormVa,&enormVr)); 286818efac9SLisandro Dalcin if (wnormtype == NORM_2) 287818efac9SLisandro Dalcin *wlte = PetscSqrtReal(PetscSqr(enormX)/2 + PetscSqr(enormV)/2); 288818efac9SLisandro Dalcin else 289818efac9SLisandro Dalcin *wlte = PetscMax(enormX,enormV); 290818efac9SLisandro Dalcin if (order) *order = 2; 291818efac9SLisandro Dalcin PetscFunctionReturn(0); 292818efac9SLisandro Dalcin } 293818efac9SLisandro Dalcin 294818efac9SLisandro Dalcin static PetscErrorCode TSRollBack_Alpha(TS ts) 295818efac9SLisandro Dalcin { 296818efac9SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 297818efac9SLisandro Dalcin 298818efac9SLisandro Dalcin PetscFunctionBegin; 2999566063dSJacob Faibussowitsch PetscCall(VecCopy(th->X0,ts->vec_sol)); 3009566063dSJacob Faibussowitsch PetscCall(VecCopy(th->V0,ts->vec_dot)); 301818efac9SLisandro Dalcin PetscFunctionReturn(0); 302818efac9SLisandro Dalcin } 303818efac9SLisandro Dalcin 304818efac9SLisandro Dalcin /* 305818efac9SLisandro Dalcin static PetscErrorCode TSInterpolate_Alpha(TS ts,PetscReal t,Vec X,Vec V) 306818efac9SLisandro Dalcin { 307818efac9SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 308818efac9SLisandro Dalcin PetscReal dt = t - ts->ptime; 309818efac9SLisandro Dalcin 310818efac9SLisandro Dalcin PetscFunctionBegin; 3119566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_dot,V)); 3129566063dSJacob Faibussowitsch PetscCall(VecAXPY(V,dt*(1-th->Gamma),th->A0)); 3139566063dSJacob Faibussowitsch PetscCall(VecAXPY(V,dt*th->Gamma,th->A1)); 3149566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol,X)); 3159566063dSJacob Faibussowitsch PetscCall(VecAXPY(X,dt,V)); 3169566063dSJacob Faibussowitsch PetscCall(VecAXPY(X,dt*dt*((PetscReal)0.5-th->Beta),th->A0)); 3179566063dSJacob Faibussowitsch PetscCall(VecAXPY(X,dt*dt*th->Beta,th->A1)); 318818efac9SLisandro Dalcin PetscFunctionReturn(0); 319818efac9SLisandro Dalcin } 320818efac9SLisandro Dalcin */ 321818efac9SLisandro Dalcin 322818efac9SLisandro Dalcin static PetscErrorCode SNESTSFormFunction_Alpha(PETSC_UNUSED SNES snes,Vec X,Vec F,TS ts) 323818efac9SLisandro Dalcin { 324818efac9SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 325818efac9SLisandro Dalcin PetscReal ta = th->stage_time; 326818efac9SLisandro Dalcin Vec Xa = th->Xa, Va = th->Va, Aa = th->Aa; 327818efac9SLisandro Dalcin 328818efac9SLisandro Dalcin PetscFunctionBegin; 3299566063dSJacob Faibussowitsch PetscCall(TSAlpha_StageVecs(ts,X)); 330818efac9SLisandro Dalcin /* F = Function(ta,Xa,Va,Aa) */ 3319566063dSJacob Faibussowitsch PetscCall(TSComputeI2Function(ts,ta,Xa,Va,Aa,F)); 3329566063dSJacob Faibussowitsch PetscCall(VecScale(F,th->scale_F)); 333818efac9SLisandro Dalcin PetscFunctionReturn(0); 334818efac9SLisandro Dalcin } 335818efac9SLisandro Dalcin 336818efac9SLisandro Dalcin static PetscErrorCode SNESTSFormJacobian_Alpha(PETSC_UNUSED SNES snes,PETSC_UNUSED Vec X,Mat J,Mat P,TS ts) 337818efac9SLisandro Dalcin { 338818efac9SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 339818efac9SLisandro Dalcin PetscReal ta = th->stage_time; 340818efac9SLisandro Dalcin Vec Xa = th->Xa, Va = th->Va, Aa = th->Aa; 341818efac9SLisandro Dalcin PetscReal dVdX = th->shift_V, dAdX = th->shift_A; 342818efac9SLisandro Dalcin 343818efac9SLisandro Dalcin PetscFunctionBegin; 344818efac9SLisandro Dalcin /* J,P = Jacobian(ta,Xa,Va,Aa) */ 3459566063dSJacob Faibussowitsch PetscCall(TSComputeI2Jacobian(ts,ta,Xa,Va,Aa,dVdX,dAdX,J,P)); 346818efac9SLisandro Dalcin PetscFunctionReturn(0); 347818efac9SLisandro Dalcin } 348818efac9SLisandro Dalcin 349818efac9SLisandro Dalcin static PetscErrorCode TSReset_Alpha(TS ts) 350818efac9SLisandro Dalcin { 351818efac9SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 352818efac9SLisandro Dalcin 353818efac9SLisandro Dalcin PetscFunctionBegin; 3549566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->X0)); 3559566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->Xa)); 3569566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->X1)); 3579566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->V0)); 3589566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->Va)); 3599566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->V1)); 3609566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->A0)); 3619566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->Aa)); 3629566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->A1)); 3639566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->vec_sol_prev)); 3649566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->vec_dot_prev)); 3659566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->vec_lte_work[0])); 3669566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->vec_lte_work[1])); 367818efac9SLisandro Dalcin PetscFunctionReturn(0); 368818efac9SLisandro Dalcin } 369818efac9SLisandro Dalcin 370818efac9SLisandro Dalcin static PetscErrorCode TSDestroy_Alpha(TS ts) 371818efac9SLisandro Dalcin { 372818efac9SLisandro Dalcin PetscFunctionBegin; 3739566063dSJacob Faibussowitsch PetscCall(TSReset_Alpha(ts)); 3749566063dSJacob Faibussowitsch PetscCall(PetscFree(ts->data)); 375818efac9SLisandro Dalcin 3769566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSAlpha2SetRadius_C",NULL)); 3779566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSAlpha2SetParams_C",NULL)); 3789566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSAlpha2GetParams_C",NULL)); 379818efac9SLisandro Dalcin PetscFunctionReturn(0); 380818efac9SLisandro Dalcin } 381818efac9SLisandro Dalcin 382818efac9SLisandro Dalcin static PetscErrorCode TSSetUp_Alpha(TS ts) 383818efac9SLisandro Dalcin { 384818efac9SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 3852ffb9264SLisandro Dalcin PetscBool match; 386818efac9SLisandro Dalcin 387818efac9SLisandro Dalcin PetscFunctionBegin; 3889566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol,&th->X0)); 3899566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol,&th->Xa)); 3909566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol,&th->X1)); 3919566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol,&th->V0)); 3929566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol,&th->Va)); 3939566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol,&th->V1)); 3949566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol,&th->A0)); 3959566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol,&th->Aa)); 3969566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol,&th->A1)); 397818efac9SLisandro Dalcin 3989566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts,&ts->adapt)); 3999566063dSJacob Faibussowitsch PetscCall(TSAdaptCandidatesClear(ts->adapt)); 4009566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)ts->adapt,TSADAPTNONE,&match)); 4012ffb9264SLisandro Dalcin if (!match) { 4029566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol,&th->vec_sol_prev)); 4039566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol,&th->vec_dot_prev)); 4049566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol,&th->vec_lte_work[0])); 4059566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol,&th->vec_lte_work[1])); 406818efac9SLisandro Dalcin } 407818efac9SLisandro Dalcin 4089566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts,&ts->snes)); 409818efac9SLisandro Dalcin PetscFunctionReturn(0); 410818efac9SLisandro Dalcin } 411818efac9SLisandro Dalcin 412818efac9SLisandro Dalcin static PetscErrorCode TSSetFromOptions_Alpha(PetscOptionItems *PetscOptionsObject,TS ts) 413818efac9SLisandro Dalcin { 414818efac9SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 415818efac9SLisandro Dalcin 416818efac9SLisandro Dalcin PetscFunctionBegin; 4179566063dSJacob Faibussowitsch PetscCall(PetscOptionsHead(PetscOptionsObject,"Generalized-Alpha ODE solver options")); 418818efac9SLisandro Dalcin { 419818efac9SLisandro Dalcin PetscBool flg; 420818efac9SLisandro Dalcin PetscReal radius = 1; 4219566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_alpha_radius","Spectral radius (high-frequency dissipation)","TSAlpha2SetRadius",radius,&radius,&flg)); 4229566063dSJacob Faibussowitsch if (flg) PetscCall(TSAlpha2SetRadius(ts,radius)); 4239566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_alpha_alpha_m","Algorithmic parameter alpha_m","TSAlpha2SetParams",th->Alpha_m,&th->Alpha_m,NULL)); 4249566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_alpha_alpha_f","Algorithmic parameter alpha_f","TSAlpha2SetParams",th->Alpha_f,&th->Alpha_f,NULL)); 4259566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_alpha_gamma","Algorithmic parameter gamma","TSAlpha2SetParams",th->Gamma,&th->Gamma,NULL)); 4269566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_alpha_beta","Algorithmic parameter beta","TSAlpha2SetParams",th->Beta,&th->Beta,NULL)); 4279566063dSJacob Faibussowitsch PetscCall(TSAlpha2SetParams(ts,th->Alpha_m,th->Alpha_f,th->Gamma,th->Beta)); 428818efac9SLisandro Dalcin } 4299566063dSJacob Faibussowitsch PetscCall(PetscOptionsTail()); 430818efac9SLisandro Dalcin PetscFunctionReturn(0); 431818efac9SLisandro Dalcin } 432818efac9SLisandro Dalcin 433818efac9SLisandro Dalcin static PetscErrorCode TSView_Alpha(TS ts,PetscViewer viewer) 434818efac9SLisandro Dalcin { 435818efac9SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 436818efac9SLisandro Dalcin PetscBool iascii; 437818efac9SLisandro Dalcin 438818efac9SLisandro Dalcin PetscFunctionBegin; 4399566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 440818efac9SLisandro Dalcin if (iascii) { 4419566063dSJacob Faibussowitsch 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)); 442818efac9SLisandro Dalcin } 443818efac9SLisandro Dalcin PetscFunctionReturn(0); 444818efac9SLisandro Dalcin } 445818efac9SLisandro Dalcin 446818efac9SLisandro Dalcin static PetscErrorCode TSAlpha2SetRadius_Alpha(TS ts,PetscReal radius) 447818efac9SLisandro Dalcin { 448818efac9SLisandro Dalcin PetscReal alpha_m,alpha_f,gamma,beta; 449818efac9SLisandro Dalcin 450818efac9SLisandro Dalcin PetscFunctionBegin; 4512c71b3e2SJacob Faibussowitsch PetscCheckFalse(radius < 0 || radius > 1,PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Radius %g not in range [0,1]",(double)radius); 452818efac9SLisandro Dalcin alpha_m = (2-radius)/(1+radius); 453818efac9SLisandro Dalcin alpha_f = 1/(1+radius); 454818efac9SLisandro Dalcin gamma = (PetscReal)0.5 + alpha_m - alpha_f; 455818efac9SLisandro Dalcin beta = (PetscReal)0.5 * (1 + alpha_m - alpha_f); beta *= beta; 4569566063dSJacob Faibussowitsch PetscCall(TSAlpha2SetParams(ts,alpha_m,alpha_f,gamma,beta)); 457818efac9SLisandro Dalcin PetscFunctionReturn(0); 458818efac9SLisandro Dalcin } 459818efac9SLisandro Dalcin 460818efac9SLisandro Dalcin static PetscErrorCode TSAlpha2SetParams_Alpha(TS ts,PetscReal alpha_m,PetscReal alpha_f,PetscReal gamma,PetscReal beta) 461818efac9SLisandro Dalcin { 462818efac9SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 463818efac9SLisandro Dalcin PetscReal tol = 100*PETSC_MACHINE_EPSILON; 464818efac9SLisandro Dalcin PetscReal res = ((PetscReal)0.5 + alpha_m - alpha_f) - gamma; 465818efac9SLisandro Dalcin 466818efac9SLisandro Dalcin PetscFunctionBegin; 467818efac9SLisandro Dalcin th->Alpha_m = alpha_m; 468818efac9SLisandro Dalcin th->Alpha_f = alpha_f; 469818efac9SLisandro Dalcin th->Gamma = gamma; 470818efac9SLisandro Dalcin th->Beta = beta; 471818efac9SLisandro Dalcin th->order = (PetscAbsReal(res) < tol) ? 2 : 1; 472818efac9SLisandro Dalcin PetscFunctionReturn(0); 473818efac9SLisandro Dalcin } 474818efac9SLisandro Dalcin 475818efac9SLisandro Dalcin static PetscErrorCode TSAlpha2GetParams_Alpha(TS ts,PetscReal *alpha_m,PetscReal *alpha_f,PetscReal *gamma,PetscReal *beta) 476818efac9SLisandro Dalcin { 477818efac9SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 478818efac9SLisandro Dalcin 479818efac9SLisandro Dalcin PetscFunctionBegin; 480818efac9SLisandro Dalcin if (alpha_m) *alpha_m = th->Alpha_m; 481818efac9SLisandro Dalcin if (alpha_f) *alpha_f = th->Alpha_f; 482818efac9SLisandro Dalcin if (gamma) *gamma = th->Gamma; 483818efac9SLisandro Dalcin if (beta) *beta = th->Beta; 484818efac9SLisandro Dalcin PetscFunctionReturn(0); 485818efac9SLisandro Dalcin } 486818efac9SLisandro Dalcin 487818efac9SLisandro Dalcin /*MC 488818efac9SLisandro Dalcin TSALPHA2 - ODE/DAE solver using the implicit Generalized-Alpha method 489818efac9SLisandro Dalcin for second-order systems 490818efac9SLisandro Dalcin 491818efac9SLisandro Dalcin Level: beginner 492818efac9SLisandro Dalcin 493818efac9SLisandro Dalcin References: 494606c0280SSatish Balay . * - J. Chung, G.M.Hubert. "A Time Integration Algorithm for Structural 495818efac9SLisandro Dalcin Dynamics with Improved Numerical Dissipation: The Generalized-alpha 496818efac9SLisandro Dalcin Method" ASME Journal of Applied Mechanics, 60, 371:375, 1993. 497818efac9SLisandro Dalcin 498818efac9SLisandro Dalcin .seealso: TS, TSCreate(), TSSetType(), TSAlpha2SetRadius(), TSAlpha2SetParams() 499818efac9SLisandro Dalcin M*/ 500818efac9SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSCreate_Alpha2(TS ts) 501818efac9SLisandro Dalcin { 502818efac9SLisandro Dalcin TS_Alpha *th; 503818efac9SLisandro Dalcin 504818efac9SLisandro Dalcin PetscFunctionBegin; 505818efac9SLisandro Dalcin ts->ops->reset = TSReset_Alpha; 506818efac9SLisandro Dalcin ts->ops->destroy = TSDestroy_Alpha; 507818efac9SLisandro Dalcin ts->ops->view = TSView_Alpha; 508818efac9SLisandro Dalcin ts->ops->setup = TSSetUp_Alpha; 509818efac9SLisandro Dalcin ts->ops->setfromoptions = TSSetFromOptions_Alpha; 510818efac9SLisandro Dalcin ts->ops->step = TSStep_Alpha; 511818efac9SLisandro Dalcin ts->ops->evaluatewlte = TSEvaluateWLTE_Alpha; 512818efac9SLisandro Dalcin ts->ops->rollback = TSRollBack_Alpha; 513818efac9SLisandro Dalcin /*ts->ops->interpolate = TSInterpolate_Alpha;*/ 514818efac9SLisandro Dalcin ts->ops->snesfunction = SNESTSFormFunction_Alpha; 515818efac9SLisandro Dalcin ts->ops->snesjacobian = SNESTSFormJacobian_Alpha; 5162ffb9264SLisandro Dalcin ts->default_adapt_type = TSADAPTNONE; 517818efac9SLisandro Dalcin 518efd4aadfSBarry Smith ts->usessnes = PETSC_TRUE; 519efd4aadfSBarry Smith 5209566063dSJacob Faibussowitsch PetscCall(PetscNewLog(ts,&th)); 521818efac9SLisandro Dalcin ts->data = (void*)th; 522818efac9SLisandro Dalcin 523818efac9SLisandro Dalcin th->Alpha_m = 0.5; 524818efac9SLisandro Dalcin th->Alpha_f = 0.5; 525818efac9SLisandro Dalcin th->Gamma = 0.5; 526818efac9SLisandro Dalcin th->Beta = 0.25; 527818efac9SLisandro Dalcin th->order = 2; 528818efac9SLisandro Dalcin 5299566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSAlpha2SetRadius_C",TSAlpha2SetRadius_Alpha)); 5309566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSAlpha2SetParams_C",TSAlpha2SetParams_Alpha)); 5319566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSAlpha2GetParams_C",TSAlpha2GetParams_Alpha)); 532818efac9SLisandro Dalcin PetscFunctionReturn(0); 533818efac9SLisandro Dalcin } 534818efac9SLisandro Dalcin 535818efac9SLisandro Dalcin /*@ 536818efac9SLisandro Dalcin TSAlpha2SetRadius - sets the desired spectral radius of the method 537818efac9SLisandro Dalcin (i.e. high-frequency numerical damping) 538818efac9SLisandro Dalcin 539818efac9SLisandro Dalcin Logically Collective on TS 540818efac9SLisandro Dalcin 541818efac9SLisandro Dalcin The algorithmic parameters \alpha_m and \alpha_f of the 542818efac9SLisandro Dalcin generalized-\alpha method can be computed in terms of a specified 543818efac9SLisandro Dalcin spectral radius \rho in [0,1] for infinite time step in order to 544818efac9SLisandro Dalcin control high-frequency numerical damping: 545818efac9SLisandro Dalcin \alpha_m = (2-\rho)/(1+\rho) 546818efac9SLisandro Dalcin \alpha_f = 1/(1+\rho) 547818efac9SLisandro Dalcin 548d8d19677SJose E. Roman Input Parameters: 549818efac9SLisandro Dalcin + ts - timestepping context 550818efac9SLisandro Dalcin - radius - the desired spectral radius 551818efac9SLisandro Dalcin 552818efac9SLisandro Dalcin Options Database: 55367b8a455SSatish Balay . -ts_alpha_radius <radius> - set the desired spectral radius 554818efac9SLisandro Dalcin 555818efac9SLisandro Dalcin Level: intermediate 556818efac9SLisandro Dalcin 557818efac9SLisandro Dalcin .seealso: TSAlpha2SetParams(), TSAlpha2GetParams() 558818efac9SLisandro Dalcin @*/ 559818efac9SLisandro Dalcin PetscErrorCode TSAlpha2SetRadius(TS ts,PetscReal radius) 560818efac9SLisandro Dalcin { 561818efac9SLisandro Dalcin PetscFunctionBegin; 562818efac9SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 563818efac9SLisandro Dalcin PetscValidLogicalCollectiveReal(ts,radius,2); 5642c71b3e2SJacob Faibussowitsch PetscCheckFalse(radius < 0 || radius > 1,((PetscObject)ts)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Radius %g not in range [0,1]",(double)radius); 565*cac4c232SBarry Smith PetscTryMethod(ts,"TSAlpha2SetRadius_C",(TS,PetscReal),(ts,radius)); 566818efac9SLisandro Dalcin PetscFunctionReturn(0); 567818efac9SLisandro Dalcin } 568818efac9SLisandro Dalcin 569818efac9SLisandro Dalcin /*@ 570818efac9SLisandro Dalcin TSAlpha2SetParams - sets the algorithmic parameters for TSALPHA2 571818efac9SLisandro Dalcin 572818efac9SLisandro Dalcin Logically Collective on TS 573818efac9SLisandro Dalcin 574818efac9SLisandro Dalcin Second-order accuracy can be obtained so long as: 575818efac9SLisandro Dalcin \gamma = 1/2 + alpha_m - alpha_f 576818efac9SLisandro Dalcin \beta = 1/4 (1 + alpha_m - alpha_f)^2 577818efac9SLisandro Dalcin 578818efac9SLisandro Dalcin Unconditional stability requires: 579818efac9SLisandro Dalcin \alpha_m >= \alpha_f >= 1/2 580818efac9SLisandro Dalcin 581d8d19677SJose E. Roman Input Parameters: 582818efac9SLisandro Dalcin + ts - timestepping context 583a5b23f4aSJose E. Roman . \alpha_m - algorithmic parameter 584a5b23f4aSJose E. Roman . \alpha_f - algorithmic parameter 585a5b23f4aSJose E. Roman . \gamma - algorithmic parameter 586a5b23f4aSJose E. Roman - \beta - algorithmic parameter 587818efac9SLisandro Dalcin 588818efac9SLisandro Dalcin Options Database: 58967b8a455SSatish Balay + -ts_alpha_alpha_m <alpha_m> - set alpha_m 59067b8a455SSatish Balay . -ts_alpha_alpha_f <alpha_f> - set alpha_f 59167b8a455SSatish Balay . -ts_alpha_gamma <gamma> - set gamma 59267b8a455SSatish Balay - -ts_alpha_beta <beta> - set beta 593818efac9SLisandro Dalcin 594818efac9SLisandro Dalcin Note: 595818efac9SLisandro Dalcin Use of this function is normally only required to hack TSALPHA2 to 596818efac9SLisandro Dalcin use a modified integration scheme. Users should call 597818efac9SLisandro Dalcin TSAlpha2SetRadius() to set the desired spectral radius of the methods 598818efac9SLisandro Dalcin (i.e. high-frequency damping) in order so select optimal values for 599818efac9SLisandro Dalcin these parameters. 600818efac9SLisandro Dalcin 601818efac9SLisandro Dalcin Level: advanced 602818efac9SLisandro Dalcin 603818efac9SLisandro Dalcin .seealso: TSAlpha2SetRadius(), TSAlpha2GetParams() 604818efac9SLisandro Dalcin @*/ 605818efac9SLisandro Dalcin PetscErrorCode TSAlpha2SetParams(TS ts,PetscReal alpha_m,PetscReal alpha_f,PetscReal gamma,PetscReal beta) 606818efac9SLisandro Dalcin { 607818efac9SLisandro Dalcin PetscFunctionBegin; 608818efac9SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 609818efac9SLisandro Dalcin PetscValidLogicalCollectiveReal(ts,alpha_m,2); 610818efac9SLisandro Dalcin PetscValidLogicalCollectiveReal(ts,alpha_f,3); 611818efac9SLisandro Dalcin PetscValidLogicalCollectiveReal(ts,gamma,4); 612818efac9SLisandro Dalcin PetscValidLogicalCollectiveReal(ts,beta,5); 613*cac4c232SBarry Smith PetscTryMethod(ts,"TSAlpha2SetParams_C",(TS,PetscReal,PetscReal,PetscReal,PetscReal),(ts,alpha_m,alpha_f,gamma,beta)); 614818efac9SLisandro Dalcin PetscFunctionReturn(0); 615818efac9SLisandro Dalcin } 616818efac9SLisandro Dalcin 617818efac9SLisandro Dalcin /*@ 618818efac9SLisandro Dalcin TSAlpha2GetParams - gets the algorithmic parameters for TSALPHA2 619818efac9SLisandro Dalcin 620818efac9SLisandro Dalcin Not Collective 621818efac9SLisandro Dalcin 622818efac9SLisandro Dalcin Input Parameter: 623818efac9SLisandro Dalcin . ts - timestepping context 624818efac9SLisandro Dalcin 625818efac9SLisandro Dalcin Output Parameters: 626818efac9SLisandro Dalcin + \alpha_m - algorithmic parameter 627818efac9SLisandro Dalcin . \alpha_f - algorithmic parameter 628818efac9SLisandro Dalcin . \gamma - algorithmic parameter 629818efac9SLisandro Dalcin - \beta - algorithmic parameter 630818efac9SLisandro Dalcin 631818efac9SLisandro Dalcin Note: 632818efac9SLisandro Dalcin Use of this function is normally only required to hack TSALPHA2 to 633818efac9SLisandro Dalcin use a modified integration scheme. Users should call 634818efac9SLisandro Dalcin TSAlpha2SetRadius() to set the high-frequency damping (i.e. spectral 635818efac9SLisandro Dalcin radius of the method) in order so select optimal values for these 636818efac9SLisandro Dalcin parameters. 637818efac9SLisandro Dalcin 638818efac9SLisandro Dalcin Level: advanced 639818efac9SLisandro Dalcin 640818efac9SLisandro Dalcin .seealso: TSAlpha2SetRadius(), TSAlpha2SetParams() 641818efac9SLisandro Dalcin @*/ 642818efac9SLisandro Dalcin PetscErrorCode TSAlpha2GetParams(TS ts,PetscReal *alpha_m,PetscReal *alpha_f,PetscReal *gamma,PetscReal *beta) 643818efac9SLisandro Dalcin { 644818efac9SLisandro Dalcin PetscFunctionBegin; 645818efac9SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 646818efac9SLisandro Dalcin if (alpha_m) PetscValidRealPointer(alpha_m,2); 647818efac9SLisandro Dalcin if (alpha_f) PetscValidRealPointer(alpha_f,3); 648818efac9SLisandro Dalcin if (gamma) PetscValidRealPointer(gamma,4); 649818efac9SLisandro Dalcin if (beta) PetscValidRealPointer(beta,5); 650*cac4c232SBarry Smith PetscUseMethod(ts,"TSAlpha2GetParams_C",(TS,PetscReal*,PetscReal*,PetscReal*,PetscReal*),(ts,alpha_m,alpha_f,gamma,beta)); 651818efac9SLisandro Dalcin PetscFunctionReturn(0); 652818efac9SLisandro Dalcin } 653