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 PetscErrorCode ierr; 74818efac9SLisandro Dalcin 75818efac9SLisandro Dalcin PetscFunctionBegin; 76818efac9SLisandro Dalcin /* A1 = ... */ 77818efac9SLisandro Dalcin ierr = VecWAXPY(A1,-1.0,X0,X1);CHKERRQ(ierr); 78818efac9SLisandro Dalcin ierr = VecAXPY (A1,-dt,V0);CHKERRQ(ierr); 79818efac9SLisandro Dalcin ierr = VecAXPBY(A1,-(1-2*Beta)/(2*Beta),1/(dt*dt*Beta),A0);CHKERRQ(ierr); 80818efac9SLisandro Dalcin /* V1 = ... */ 81818efac9SLisandro Dalcin ierr = VecWAXPY(V1,(1.0-Gamma)/Gamma,A0,A1);CHKERRQ(ierr); 82818efac9SLisandro Dalcin ierr = VecAYPX (V1,dt*Gamma,V0);CHKERRQ(ierr); 83818efac9SLisandro Dalcin /* Xa = X0 + Alpha_f*(X1-X0) */ 84818efac9SLisandro Dalcin ierr = VecWAXPY(Xa,-1.0,X0,X1);CHKERRQ(ierr); 85818efac9SLisandro Dalcin ierr = VecAYPX (Xa,Alpha_f,X0);CHKERRQ(ierr); 86818efac9SLisandro Dalcin /* Va = V0 + Alpha_f*(V1-V0) */ 87818efac9SLisandro Dalcin ierr = VecWAXPY(Va,-1.0,V0,V1);CHKERRQ(ierr); 88818efac9SLisandro Dalcin ierr = VecAYPX (Va,Alpha_f,V0);CHKERRQ(ierr); 89818efac9SLisandro Dalcin /* Aa = A0 + Alpha_m*(A1-A0) */ 90818efac9SLisandro Dalcin ierr = VecWAXPY(Aa,-1.0,A0,A1);CHKERRQ(ierr); 91818efac9SLisandro Dalcin ierr = VecAYPX (Aa,Alpha_m,A0);CHKERRQ(ierr); 92818efac9SLisandro Dalcin PetscFunctionReturn(0); 93818efac9SLisandro Dalcin } 94818efac9SLisandro Dalcin 95818efac9SLisandro Dalcin static PetscErrorCode TS_SNESSolve(TS ts,Vec b,Vec x) 96818efac9SLisandro Dalcin { 97818efac9SLisandro Dalcin PetscInt nits,lits; 98818efac9SLisandro Dalcin PetscErrorCode ierr; 99818efac9SLisandro Dalcin 100818efac9SLisandro Dalcin PetscFunctionBegin; 101818efac9SLisandro Dalcin ierr = SNESSolve(ts->snes,b,x);CHKERRQ(ierr); 102818efac9SLisandro Dalcin ierr = SNESGetIterationNumber(ts->snes,&nits);CHKERRQ(ierr); 103818efac9SLisandro Dalcin ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr); 104818efac9SLisandro Dalcin ts->snes_its += nits; ts->ksp_its += lits; 105818efac9SLisandro Dalcin PetscFunctionReturn(0); 106818efac9SLisandro Dalcin } 107818efac9SLisandro Dalcin 108818efac9SLisandro Dalcin /* 109818efac9SLisandro Dalcin Compute a consistent initial state for the generalized-alpha method. 110818efac9SLisandro Dalcin - Solve two successive backward Euler steps with halved time step. 111818efac9SLisandro Dalcin - Compute the initial second time derivative using backward differences. 112818efac9SLisandro Dalcin - If using adaptivity, estimate the LTE of the initial step. 113818efac9SLisandro Dalcin */ 114818efac9SLisandro Dalcin static PetscErrorCode TSAlpha_Restart(TS ts,PetscBool *initok) 115818efac9SLisandro Dalcin { 116818efac9SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 117818efac9SLisandro Dalcin PetscReal time_step; 118818efac9SLisandro Dalcin PetscReal alpha_m,alpha_f,gamma,beta; 119818efac9SLisandro Dalcin Vec X0 = ts->vec_sol, X1, X2 = th->X1; 120818efac9SLisandro Dalcin Vec V0 = ts->vec_dot, V1, V2 = th->V1; 121818efac9SLisandro Dalcin PetscBool stageok; 122818efac9SLisandro Dalcin PetscErrorCode ierr; 123818efac9SLisandro Dalcin 124818efac9SLisandro Dalcin PetscFunctionBegin; 125818efac9SLisandro Dalcin ierr = VecDuplicate(X0,&X1);CHKERRQ(ierr); 126818efac9SLisandro Dalcin ierr = VecDuplicate(V0,&V1);CHKERRQ(ierr); 127818efac9SLisandro Dalcin 128818efac9SLisandro Dalcin /* Setup backward Euler with halved time step */ 129818efac9SLisandro Dalcin ierr = TSAlpha2GetParams(ts,&alpha_m,&alpha_f,&gamma,&beta);CHKERRQ(ierr); 130818efac9SLisandro Dalcin ierr = TSAlpha2SetParams(ts,1,1,1,0.5);CHKERRQ(ierr); 131818efac9SLisandro Dalcin ierr = TSGetTimeStep(ts,&time_step);CHKERRQ(ierr); 132818efac9SLisandro Dalcin ts->time_step = time_step/2; 133818efac9SLisandro Dalcin ierr = TSAlpha_StageTime(ts);CHKERRQ(ierr); 134818efac9SLisandro Dalcin th->stage_time = ts->ptime; 135818efac9SLisandro Dalcin ierr = VecZeroEntries(th->A0);CHKERRQ(ierr); 136818efac9SLisandro Dalcin 137818efac9SLisandro Dalcin /* First BE step, (t0,X0,V0) -> (t1,X1,V1) */ 138818efac9SLisandro Dalcin th->stage_time += ts->time_step; 139818efac9SLisandro Dalcin ierr = VecCopy(X0,th->X0);CHKERRQ(ierr); 140818efac9SLisandro Dalcin ierr = VecCopy(V0,th->V0);CHKERRQ(ierr); 141818efac9SLisandro Dalcin ierr = TSPreStage(ts,th->stage_time);CHKERRQ(ierr); 142818efac9SLisandro Dalcin ierr = VecCopy(th->X0,X1);CHKERRQ(ierr); 143818efac9SLisandro Dalcin ierr = TS_SNESSolve(ts,NULL,X1);CHKERRQ(ierr); 144818efac9SLisandro Dalcin ierr = VecCopy(th->V1,V1);CHKERRQ(ierr); 145818efac9SLisandro Dalcin ierr = TSPostStage(ts,th->stage_time,0,&X1);CHKERRQ(ierr); 146818efac9SLisandro Dalcin ierr = TSAdaptCheckStage(ts->adapt,ts,th->stage_time,X1,&stageok);CHKERRQ(ierr); 147818efac9SLisandro Dalcin if (!stageok) goto finally; 148818efac9SLisandro Dalcin 149818efac9SLisandro Dalcin /* Second BE step, (t1,X1,V1) -> (t2,X2,V2) */ 150818efac9SLisandro Dalcin th->stage_time += ts->time_step; 151818efac9SLisandro Dalcin ierr = VecCopy(X1,th->X0);CHKERRQ(ierr); 152818efac9SLisandro Dalcin ierr = VecCopy(V1,th->V0);CHKERRQ(ierr); 153818efac9SLisandro Dalcin ierr = TSPreStage(ts,th->stage_time);CHKERRQ(ierr); 154818efac9SLisandro Dalcin ierr = VecCopy(th->X0,X2);CHKERRQ(ierr); 155818efac9SLisandro Dalcin ierr = TS_SNESSolve(ts,NULL,X2);CHKERRQ(ierr); 156818efac9SLisandro Dalcin ierr = VecCopy(th->V1,V2);CHKERRQ(ierr); 157818efac9SLisandro Dalcin ierr = TSPostStage(ts,th->stage_time,0,&X2);CHKERRQ(ierr); 158818efac9SLisandro Dalcin ierr = TSAdaptCheckStage(ts->adapt,ts,th->stage_time,X1,&stageok);CHKERRQ(ierr); 159818efac9SLisandro Dalcin if (!stageok) goto finally; 160818efac9SLisandro Dalcin 161818efac9SLisandro Dalcin /* Compute A0 ~ dV/dt at t0 with backward differences */ 162818efac9SLisandro Dalcin ierr = VecZeroEntries(th->A0);CHKERRQ(ierr); 163818efac9SLisandro Dalcin ierr = VecAXPY(th->A0,-3/ts->time_step,V0);CHKERRQ(ierr); 164818efac9SLisandro Dalcin ierr = VecAXPY(th->A0,+4/ts->time_step,V1);CHKERRQ(ierr); 165818efac9SLisandro Dalcin ierr = VecAXPY(th->A0,-1/ts->time_step,V2);CHKERRQ(ierr); 166818efac9SLisandro Dalcin 167818efac9SLisandro Dalcin /* Rough, lower-order estimate LTE of the initial step */ 1682ffb9264SLisandro Dalcin if (th->vec_lte_work[0]) { 169818efac9SLisandro Dalcin ierr = VecZeroEntries(th->vec_lte_work[0]);CHKERRQ(ierr); 170818efac9SLisandro Dalcin ierr = VecAXPY(th->vec_lte_work[0],+2,X2);CHKERRQ(ierr); 171818efac9SLisandro Dalcin ierr = VecAXPY(th->vec_lte_work[0],-4,X1);CHKERRQ(ierr); 172818efac9SLisandro Dalcin ierr = VecAXPY(th->vec_lte_work[0],+2,X0);CHKERRQ(ierr); 173818efac9SLisandro Dalcin } 1742ffb9264SLisandro Dalcin if (th->vec_lte_work[1]) { 175818efac9SLisandro Dalcin ierr = VecZeroEntries(th->vec_lte_work[1]);CHKERRQ(ierr); 176818efac9SLisandro Dalcin ierr = VecAXPY(th->vec_lte_work[1],+2,V2);CHKERRQ(ierr); 177818efac9SLisandro Dalcin ierr = VecAXPY(th->vec_lte_work[1],-4,V1);CHKERRQ(ierr); 178818efac9SLisandro Dalcin ierr = VecAXPY(th->vec_lte_work[1],+2,V0);CHKERRQ(ierr); 179818efac9SLisandro Dalcin } 180818efac9SLisandro Dalcin 181818efac9SLisandro Dalcin finally: 182818efac9SLisandro Dalcin /* Revert TSAlpha to the initial state (t0,X0,V0) */ 183818efac9SLisandro Dalcin if (initok) *initok = stageok; 184818efac9SLisandro Dalcin ierr = TSSetTimeStep(ts,time_step);CHKERRQ(ierr); 185818efac9SLisandro Dalcin ierr = TSAlpha2SetParams(ts,alpha_m,alpha_f,gamma,beta);CHKERRQ(ierr); 186818efac9SLisandro Dalcin ierr = VecCopy(ts->vec_sol,th->X0);CHKERRQ(ierr); 187818efac9SLisandro Dalcin ierr = VecCopy(ts->vec_dot,th->V0);CHKERRQ(ierr); 188818efac9SLisandro Dalcin 189818efac9SLisandro Dalcin ierr = VecDestroy(&X1);CHKERRQ(ierr); 190818efac9SLisandro Dalcin ierr = VecDestroy(&V1);CHKERRQ(ierr); 191818efac9SLisandro Dalcin PetscFunctionReturn(0); 192818efac9SLisandro Dalcin } 193818efac9SLisandro Dalcin 194818efac9SLisandro Dalcin static PetscErrorCode TSStep_Alpha(TS ts) 195818efac9SLisandro Dalcin { 196818efac9SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 197818efac9SLisandro Dalcin PetscInt rejections = 0; 198818efac9SLisandro Dalcin PetscBool stageok,accept = PETSC_TRUE; 199818efac9SLisandro Dalcin PetscReal next_time_step = ts->time_step; 200818efac9SLisandro Dalcin PetscErrorCode ierr; 201818efac9SLisandro Dalcin 202818efac9SLisandro Dalcin PetscFunctionBegin; 203818efac9SLisandro Dalcin ierr = PetscCitationsRegister(citation,&cited);CHKERRQ(ierr); 204818efac9SLisandro Dalcin 205818efac9SLisandro Dalcin if (!ts->steprollback) { 2062ffb9264SLisandro Dalcin if (th->vec_sol_prev) { ierr = VecCopy(th->X0,th->vec_sol_prev);CHKERRQ(ierr); } 2072ffb9264SLisandro Dalcin if (th->vec_dot_prev) { ierr = VecCopy(th->V0,th->vec_dot_prev);CHKERRQ(ierr); } 208818efac9SLisandro Dalcin ierr = VecCopy(ts->vec_sol,th->X0);CHKERRQ(ierr); 209818efac9SLisandro Dalcin ierr = VecCopy(ts->vec_dot,th->V0);CHKERRQ(ierr); 210818efac9SLisandro Dalcin ierr = VecCopy(th->A1,th->A0);CHKERRQ(ierr); 211818efac9SLisandro Dalcin } 212818efac9SLisandro Dalcin 213818efac9SLisandro Dalcin th->status = TS_STEP_INCOMPLETE; 214818efac9SLisandro Dalcin while (!ts->reason && th->status != TS_STEP_COMPLETE) { 215818efac9SLisandro Dalcin 216818efac9SLisandro Dalcin if (ts->steprestart) { 217818efac9SLisandro Dalcin ierr = TSAlpha_Restart(ts,&stageok);CHKERRQ(ierr); 218818efac9SLisandro Dalcin if (!stageok) goto reject_step; 219818efac9SLisandro Dalcin } 220818efac9SLisandro Dalcin 221818efac9SLisandro Dalcin ierr = TSAlpha_StageTime(ts);CHKERRQ(ierr); 222818efac9SLisandro Dalcin ierr = VecCopy(th->X0,th->X1);CHKERRQ(ierr); 223818efac9SLisandro Dalcin ierr = TSPreStage(ts,th->stage_time);CHKERRQ(ierr); 224818efac9SLisandro Dalcin ierr = TS_SNESSolve(ts,NULL,th->X1);CHKERRQ(ierr); 225818efac9SLisandro Dalcin ierr = TSPostStage(ts,th->stage_time,0,&th->Xa);CHKERRQ(ierr); 226818efac9SLisandro Dalcin ierr = TSAdaptCheckStage(ts->adapt,ts,th->stage_time,th->Xa,&stageok);CHKERRQ(ierr); 227818efac9SLisandro Dalcin if (!stageok) goto reject_step; 228818efac9SLisandro Dalcin 229818efac9SLisandro Dalcin th->status = TS_STEP_PENDING; 230818efac9SLisandro Dalcin ierr = VecCopy(th->X1,ts->vec_sol);CHKERRQ(ierr); 231818efac9SLisandro Dalcin ierr = VecCopy(th->V1,ts->vec_dot);CHKERRQ(ierr); 232818efac9SLisandro Dalcin ierr = TSAdaptChoose(ts->adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr); 233818efac9SLisandro Dalcin th->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 234818efac9SLisandro Dalcin if (!accept) { 235818efac9SLisandro Dalcin ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr); 236818efac9SLisandro Dalcin ierr = VecCopy(th->V0,ts->vec_dot);CHKERRQ(ierr); 237818efac9SLisandro Dalcin ts->time_step = next_time_step; 238818efac9SLisandro Dalcin goto reject_step; 239818efac9SLisandro Dalcin } 240818efac9SLisandro Dalcin 241818efac9SLisandro Dalcin ts->ptime += ts->time_step; 242818efac9SLisandro Dalcin ts->time_step = next_time_step; 243818efac9SLisandro Dalcin break; 244818efac9SLisandro Dalcin 245818efac9SLisandro Dalcin reject_step: 246818efac9SLisandro Dalcin ts->reject++; accept = PETSC_FALSE; 247818efac9SLisandro Dalcin if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 248818efac9SLisandro Dalcin ts->reason = TS_DIVERGED_STEP_REJECTED; 249818efac9SLisandro Dalcin ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr); 250818efac9SLisandro Dalcin } 251818efac9SLisandro Dalcin 252818efac9SLisandro Dalcin } 253818efac9SLisandro Dalcin PetscFunctionReturn(0); 254818efac9SLisandro Dalcin } 255818efac9SLisandro Dalcin 256818efac9SLisandro Dalcin static PetscErrorCode TSEvaluateWLTE_Alpha(TS ts,NormType wnormtype,PetscInt *order,PetscReal *wlte) 257818efac9SLisandro Dalcin { 258818efac9SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 259818efac9SLisandro Dalcin Vec X = th->X1; /* X = solution */ 260818efac9SLisandro Dalcin Vec V = th->V1; /* V = solution */ 261818efac9SLisandro Dalcin Vec Y = th->vec_lte_work[0]; /* Y = X + LTE */ 262818efac9SLisandro Dalcin Vec Z = th->vec_lte_work[1]; /* Z = V + LTE */ 2637453f775SEmil Constantinescu PetscReal enormX,enormV,enormXa,enormVa,enormXr,enormVr; 264818efac9SLisandro Dalcin PetscErrorCode ierr; 265818efac9SLisandro Dalcin 266818efac9SLisandro Dalcin PetscFunctionBegin; 2672ffb9264SLisandro Dalcin if (!th->vec_sol_prev) {*wlte = -1; PetscFunctionReturn(0);} 2682ffb9264SLisandro Dalcin if (!th->vec_dot_prev) {*wlte = -1; PetscFunctionReturn(0);} 2692ffb9264SLisandro Dalcin if (!th->vec_lte_work[0]) {*wlte = -1; PetscFunctionReturn(0);} 2702ffb9264SLisandro Dalcin if (!th->vec_lte_work[1]) {*wlte = -1; PetscFunctionReturn(0);} 271818efac9SLisandro Dalcin if (ts->steprestart) { 2722ffb9264SLisandro Dalcin /* th->vec_lte_prev is set to the LTE in TSAlpha_Restart() */ 273818efac9SLisandro Dalcin ierr = VecAXPY(Y,1,X);CHKERRQ(ierr); 274818efac9SLisandro Dalcin ierr = VecAXPY(Z,1,V);CHKERRQ(ierr); 275818efac9SLisandro Dalcin } else { 276818efac9SLisandro Dalcin /* Compute LTE using backward differences with non-constant time step */ 277818efac9SLisandro Dalcin PetscReal h = ts->time_step, h_prev = ts->ptime - ts->ptime_prev; 278818efac9SLisandro Dalcin PetscReal a = 1 + h_prev/h; 279818efac9SLisandro Dalcin PetscScalar scal[3]; Vec vecX[3],vecV[3]; 280818efac9SLisandro Dalcin scal[0] = +1/a; scal[1] = -1/(a-1); scal[2] = +1/(a*(a-1)); 281818efac9SLisandro Dalcin vecX[0] = th->X1; vecX[1] = th->X0; vecX[2] = th->vec_sol_prev; 282818efac9SLisandro Dalcin vecV[0] = th->V1; vecV[1] = th->V0; vecV[2] = th->vec_dot_prev; 283818efac9SLisandro Dalcin ierr = VecCopy(X,Y);CHKERRQ(ierr); 284818efac9SLisandro Dalcin ierr = VecMAXPY(Y,3,scal,vecX);CHKERRQ(ierr); 285818efac9SLisandro Dalcin ierr = VecCopy(V,Z);CHKERRQ(ierr); 286818efac9SLisandro Dalcin ierr = VecMAXPY(Z,3,scal,vecV);CHKERRQ(ierr); 287818efac9SLisandro Dalcin } 288818efac9SLisandro Dalcin /* XXX ts->atol and ts->vatol are not appropriate for computing enormV */ 2897453f775SEmil Constantinescu ierr = TSErrorWeightedNorm(ts,X,Y,wnormtype,&enormX,&enormXa,&enormXr);CHKERRQ(ierr); 2907453f775SEmil Constantinescu ierr = TSErrorWeightedNorm(ts,V,Z,wnormtype,&enormV,&enormVa,&enormVr);CHKERRQ(ierr); 291818efac9SLisandro Dalcin if (wnormtype == NORM_2) 292818efac9SLisandro Dalcin *wlte = PetscSqrtReal(PetscSqr(enormX)/2 + PetscSqr(enormV)/2); 293818efac9SLisandro Dalcin else 294818efac9SLisandro Dalcin *wlte = PetscMax(enormX,enormV); 295818efac9SLisandro Dalcin if (order) *order = 2; 296818efac9SLisandro Dalcin PetscFunctionReturn(0); 297818efac9SLisandro Dalcin } 298818efac9SLisandro Dalcin 299818efac9SLisandro Dalcin static PetscErrorCode TSRollBack_Alpha(TS ts) 300818efac9SLisandro Dalcin { 301818efac9SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 302818efac9SLisandro Dalcin PetscErrorCode ierr; 303818efac9SLisandro Dalcin 304818efac9SLisandro Dalcin PetscFunctionBegin; 305818efac9SLisandro Dalcin ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr); 306818efac9SLisandro Dalcin ierr = VecCopy(th->V0,ts->vec_dot);CHKERRQ(ierr); 307818efac9SLisandro Dalcin PetscFunctionReturn(0); 308818efac9SLisandro Dalcin } 309818efac9SLisandro Dalcin 310818efac9SLisandro Dalcin /* 311818efac9SLisandro Dalcin static PetscErrorCode TSInterpolate_Alpha(TS ts,PetscReal t,Vec X,Vec V) 312818efac9SLisandro Dalcin { 313818efac9SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 314818efac9SLisandro Dalcin PetscReal dt = t - ts->ptime; 315818efac9SLisandro Dalcin PetscErrorCode ierr; 316818efac9SLisandro Dalcin 317818efac9SLisandro Dalcin PetscFunctionBegin; 318818efac9SLisandro Dalcin ierr = VecCopy(ts->vec_dot,V);CHKERRQ(ierr); 319818efac9SLisandro Dalcin ierr = VecAXPY(V,dt*(1-th->Gamma),th->A0);CHKERRQ(ierr); 320818efac9SLisandro Dalcin ierr = VecAXPY(V,dt*th->Gamma,th->A1);CHKERRQ(ierr); 321818efac9SLisandro Dalcin ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 322818efac9SLisandro Dalcin ierr = VecAXPY(X,dt,V);CHKERRQ(ierr); 323818efac9SLisandro Dalcin ierr = VecAXPY(X,dt*dt*((PetscReal)0.5-th->Beta),th->A0);CHKERRQ(ierr); 324818efac9SLisandro Dalcin ierr = VecAXPY(X,dt*dt*th->Beta,th->A1);CHKERRQ(ierr); 325818efac9SLisandro Dalcin PetscFunctionReturn(0); 326818efac9SLisandro Dalcin } 327818efac9SLisandro Dalcin */ 328818efac9SLisandro Dalcin 329818efac9SLisandro Dalcin static PetscErrorCode SNESTSFormFunction_Alpha(PETSC_UNUSED SNES snes,Vec X,Vec F,TS ts) 330818efac9SLisandro Dalcin { 331818efac9SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 332818efac9SLisandro Dalcin PetscReal ta = th->stage_time; 333818efac9SLisandro Dalcin Vec Xa = th->Xa, Va = th->Va, Aa = th->Aa; 334818efac9SLisandro Dalcin PetscErrorCode ierr; 335818efac9SLisandro Dalcin 336818efac9SLisandro Dalcin PetscFunctionBegin; 337818efac9SLisandro Dalcin ierr = TSAlpha_StageVecs(ts,X);CHKERRQ(ierr); 338818efac9SLisandro Dalcin /* F = Function(ta,Xa,Va,Aa) */ 339818efac9SLisandro Dalcin ierr = TSComputeI2Function(ts,ta,Xa,Va,Aa,F);CHKERRQ(ierr); 340818efac9SLisandro Dalcin ierr = VecScale(F,th->scale_F);CHKERRQ(ierr); 341818efac9SLisandro Dalcin PetscFunctionReturn(0); 342818efac9SLisandro Dalcin } 343818efac9SLisandro Dalcin 344818efac9SLisandro Dalcin static PetscErrorCode SNESTSFormJacobian_Alpha(PETSC_UNUSED SNES snes,PETSC_UNUSED Vec X,Mat J,Mat P,TS ts) 345818efac9SLisandro Dalcin { 346818efac9SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 347818efac9SLisandro Dalcin PetscReal ta = th->stage_time; 348818efac9SLisandro Dalcin Vec Xa = th->Xa, Va = th->Va, Aa = th->Aa; 349818efac9SLisandro Dalcin PetscReal dVdX = th->shift_V, dAdX = th->shift_A; 350818efac9SLisandro Dalcin PetscErrorCode ierr; 351818efac9SLisandro Dalcin 352818efac9SLisandro Dalcin PetscFunctionBegin; 353818efac9SLisandro Dalcin /* J,P = Jacobian(ta,Xa,Va,Aa) */ 354818efac9SLisandro Dalcin ierr = TSComputeI2Jacobian(ts,ta,Xa,Va,Aa,dVdX,dAdX,J,P);CHKERRQ(ierr); 355818efac9SLisandro Dalcin PetscFunctionReturn(0); 356818efac9SLisandro Dalcin } 357818efac9SLisandro Dalcin 358818efac9SLisandro Dalcin static PetscErrorCode TSReset_Alpha(TS ts) 359818efac9SLisandro Dalcin { 360818efac9SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 361818efac9SLisandro Dalcin PetscErrorCode ierr; 362818efac9SLisandro Dalcin 363818efac9SLisandro Dalcin PetscFunctionBegin; 364818efac9SLisandro Dalcin ierr = VecDestroy(&th->X0);CHKERRQ(ierr); 365818efac9SLisandro Dalcin ierr = VecDestroy(&th->Xa);CHKERRQ(ierr); 366818efac9SLisandro Dalcin ierr = VecDestroy(&th->X1);CHKERRQ(ierr); 367818efac9SLisandro Dalcin ierr = VecDestroy(&th->V0);CHKERRQ(ierr); 368818efac9SLisandro Dalcin ierr = VecDestroy(&th->Va);CHKERRQ(ierr); 369818efac9SLisandro Dalcin ierr = VecDestroy(&th->V1);CHKERRQ(ierr); 370818efac9SLisandro Dalcin ierr = VecDestroy(&th->A0);CHKERRQ(ierr); 371818efac9SLisandro Dalcin ierr = VecDestroy(&th->Aa);CHKERRQ(ierr); 372818efac9SLisandro Dalcin ierr = VecDestroy(&th->A1);CHKERRQ(ierr); 373818efac9SLisandro Dalcin ierr = VecDestroy(&th->vec_sol_prev);CHKERRQ(ierr); 374818efac9SLisandro Dalcin ierr = VecDestroy(&th->vec_dot_prev);CHKERRQ(ierr); 375818efac9SLisandro Dalcin ierr = VecDestroy(&th->vec_lte_work[0]);CHKERRQ(ierr); 376818efac9SLisandro Dalcin ierr = VecDestroy(&th->vec_lte_work[1]);CHKERRQ(ierr); 377818efac9SLisandro Dalcin PetscFunctionReturn(0); 378818efac9SLisandro Dalcin } 379818efac9SLisandro Dalcin 380818efac9SLisandro Dalcin static PetscErrorCode TSDestroy_Alpha(TS ts) 381818efac9SLisandro Dalcin { 382818efac9SLisandro Dalcin PetscErrorCode ierr; 383818efac9SLisandro Dalcin 384818efac9SLisandro Dalcin PetscFunctionBegin; 385818efac9SLisandro Dalcin ierr = TSReset_Alpha(ts);CHKERRQ(ierr); 386818efac9SLisandro Dalcin ierr = PetscFree(ts->data);CHKERRQ(ierr); 387818efac9SLisandro Dalcin 388818efac9SLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)ts,"TSAlpha2SetRadius_C",NULL);CHKERRQ(ierr); 389818efac9SLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)ts,"TSAlpha2SetParams_C",NULL);CHKERRQ(ierr); 390818efac9SLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)ts,"TSAlpha2GetParams_C",NULL);CHKERRQ(ierr); 391818efac9SLisandro Dalcin PetscFunctionReturn(0); 392818efac9SLisandro Dalcin } 393818efac9SLisandro Dalcin 394818efac9SLisandro Dalcin static PetscErrorCode TSSetUp_Alpha(TS ts) 395818efac9SLisandro Dalcin { 396818efac9SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 3972ffb9264SLisandro Dalcin PetscBool match; 398818efac9SLisandro Dalcin PetscErrorCode ierr; 399818efac9SLisandro Dalcin 400818efac9SLisandro Dalcin PetscFunctionBegin; 401818efac9SLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->X0);CHKERRQ(ierr); 402818efac9SLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->Xa);CHKERRQ(ierr); 403818efac9SLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->X1);CHKERRQ(ierr); 404818efac9SLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->V0);CHKERRQ(ierr); 405818efac9SLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->Va);CHKERRQ(ierr); 406818efac9SLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->V1);CHKERRQ(ierr); 407818efac9SLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->A0);CHKERRQ(ierr); 408818efac9SLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->Aa);CHKERRQ(ierr); 409818efac9SLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->A1);CHKERRQ(ierr); 410818efac9SLisandro Dalcin 411818efac9SLisandro Dalcin ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr); 412818efac9SLisandro Dalcin ierr = TSAdaptCandidatesClear(ts->adapt);CHKERRQ(ierr); 4132ffb9264SLisandro Dalcin ierr = PetscObjectTypeCompare((PetscObject)ts->adapt,TSADAPTNONE,&match);CHKERRQ(ierr); 4142ffb9264SLisandro Dalcin if (!match) { 415818efac9SLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->vec_sol_prev);CHKERRQ(ierr); 416818efac9SLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->vec_dot_prev);CHKERRQ(ierr); 417818efac9SLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->vec_lte_work[0]);CHKERRQ(ierr); 418818efac9SLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->vec_lte_work[1]);CHKERRQ(ierr); 419818efac9SLisandro Dalcin } 420818efac9SLisandro Dalcin 421818efac9SLisandro Dalcin ierr = TSGetSNES(ts,&ts->snes);CHKERRQ(ierr); 422818efac9SLisandro Dalcin PetscFunctionReturn(0); 423818efac9SLisandro Dalcin } 424818efac9SLisandro Dalcin 425818efac9SLisandro Dalcin static PetscErrorCode TSSetFromOptions_Alpha(PetscOptionItems *PetscOptionsObject,TS ts) 426818efac9SLisandro Dalcin { 427818efac9SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 428818efac9SLisandro Dalcin PetscErrorCode ierr; 429818efac9SLisandro Dalcin 430818efac9SLisandro Dalcin PetscFunctionBegin; 431818efac9SLisandro Dalcin ierr = PetscOptionsHead(PetscOptionsObject,"Generalized-Alpha ODE solver options");CHKERRQ(ierr); 432818efac9SLisandro Dalcin { 433818efac9SLisandro Dalcin PetscBool flg; 434818efac9SLisandro Dalcin PetscReal radius = 1; 435818efac9SLisandro Dalcin ierr = PetscOptionsReal("-ts_alpha_radius","Spectral radius (high-frequency dissipation)","TSAlpha2SetRadius",radius,&radius,&flg);CHKERRQ(ierr); 436818efac9SLisandro Dalcin if (flg) {ierr = TSAlpha2SetRadius(ts,radius);CHKERRQ(ierr);} 437818efac9SLisandro Dalcin ierr = PetscOptionsReal("-ts_alpha_alpha_m","Algoritmic parameter alpha_m","TSAlpha2SetParams",th->Alpha_m,&th->Alpha_m,NULL);CHKERRQ(ierr); 438818efac9SLisandro Dalcin ierr = PetscOptionsReal("-ts_alpha_alpha_f","Algoritmic parameter alpha_f","TSAlpha2SetParams",th->Alpha_f,&th->Alpha_f,NULL);CHKERRQ(ierr); 439818efac9SLisandro Dalcin ierr = PetscOptionsReal("-ts_alpha_gamma","Algoritmic parameter gamma","TSAlpha2SetParams",th->Gamma,&th->Gamma,NULL);CHKERRQ(ierr); 440818efac9SLisandro Dalcin ierr = PetscOptionsReal("-ts_alpha_beta","Algoritmic parameter beta","TSAlpha2SetParams",th->Beta,&th->Beta,NULL);CHKERRQ(ierr); 441818efac9SLisandro Dalcin ierr = TSAlpha2SetParams(ts,th->Alpha_m,th->Alpha_f,th->Gamma,th->Beta);CHKERRQ(ierr); 442818efac9SLisandro Dalcin } 443818efac9SLisandro Dalcin ierr = PetscOptionsTail();CHKERRQ(ierr); 444818efac9SLisandro Dalcin PetscFunctionReturn(0); 445818efac9SLisandro Dalcin } 446818efac9SLisandro Dalcin 447818efac9SLisandro Dalcin static PetscErrorCode TSView_Alpha(TS ts,PetscViewer viewer) 448818efac9SLisandro Dalcin { 449818efac9SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 450818efac9SLisandro Dalcin PetscBool iascii; 451818efac9SLisandro Dalcin PetscErrorCode ierr; 452818efac9SLisandro Dalcin 453818efac9SLisandro Dalcin PetscFunctionBegin; 454818efac9SLisandro Dalcin ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 455818efac9SLisandro Dalcin if (iascii) { 456818efac9SLisandro Dalcin ierr = 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);CHKERRQ(ierr); 457818efac9SLisandro Dalcin } 458818efac9SLisandro Dalcin PetscFunctionReturn(0); 459818efac9SLisandro Dalcin } 460818efac9SLisandro Dalcin 461818efac9SLisandro Dalcin static PetscErrorCode TSAlpha2SetRadius_Alpha(TS ts,PetscReal radius) 462818efac9SLisandro Dalcin { 463818efac9SLisandro Dalcin PetscReal alpha_m,alpha_f,gamma,beta; 464818efac9SLisandro Dalcin PetscErrorCode ierr; 465818efac9SLisandro Dalcin 466818efac9SLisandro Dalcin PetscFunctionBegin; 467818efac9SLisandro Dalcin if (radius < 0 || radius > 1) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Radius %g not in range [0,1]",(double)radius); 468818efac9SLisandro Dalcin alpha_m = (2-radius)/(1+radius); 469818efac9SLisandro Dalcin alpha_f = 1/(1+radius); 470818efac9SLisandro Dalcin gamma = (PetscReal)0.5 + alpha_m - alpha_f; 471818efac9SLisandro Dalcin beta = (PetscReal)0.5 * (1 + alpha_m - alpha_f); beta *= beta; 472818efac9SLisandro Dalcin ierr = TSAlpha2SetParams(ts,alpha_m,alpha_f,gamma,beta);CHKERRQ(ierr); 473818efac9SLisandro Dalcin PetscFunctionReturn(0); 474818efac9SLisandro Dalcin } 475818efac9SLisandro Dalcin 476818efac9SLisandro Dalcin static PetscErrorCode TSAlpha2SetParams_Alpha(TS ts,PetscReal alpha_m,PetscReal alpha_f,PetscReal gamma,PetscReal beta) 477818efac9SLisandro Dalcin { 478818efac9SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 479818efac9SLisandro Dalcin PetscReal tol = 100*PETSC_MACHINE_EPSILON; 480818efac9SLisandro Dalcin PetscReal res = ((PetscReal)0.5 + alpha_m - alpha_f) - gamma; 481818efac9SLisandro Dalcin 482818efac9SLisandro Dalcin PetscFunctionBegin; 483818efac9SLisandro Dalcin th->Alpha_m = alpha_m; 484818efac9SLisandro Dalcin th->Alpha_f = alpha_f; 485818efac9SLisandro Dalcin th->Gamma = gamma; 486818efac9SLisandro Dalcin th->Beta = beta; 487818efac9SLisandro Dalcin th->order = (PetscAbsReal(res) < tol) ? 2 : 1; 488818efac9SLisandro Dalcin PetscFunctionReturn(0); 489818efac9SLisandro Dalcin } 490818efac9SLisandro Dalcin 491818efac9SLisandro Dalcin static PetscErrorCode TSAlpha2GetParams_Alpha(TS ts,PetscReal *alpha_m,PetscReal *alpha_f,PetscReal *gamma,PetscReal *beta) 492818efac9SLisandro Dalcin { 493818efac9SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 494818efac9SLisandro Dalcin 495818efac9SLisandro Dalcin PetscFunctionBegin; 496818efac9SLisandro Dalcin if (alpha_m) *alpha_m = th->Alpha_m; 497818efac9SLisandro Dalcin if (alpha_f) *alpha_f = th->Alpha_f; 498818efac9SLisandro Dalcin if (gamma) *gamma = th->Gamma; 499818efac9SLisandro Dalcin if (beta) *beta = th->Beta; 500818efac9SLisandro Dalcin PetscFunctionReturn(0); 501818efac9SLisandro Dalcin } 502818efac9SLisandro Dalcin 503818efac9SLisandro Dalcin /*MC 504818efac9SLisandro Dalcin TSALPHA2 - ODE/DAE solver using the implicit Generalized-Alpha method 505818efac9SLisandro Dalcin for second-order systems 506818efac9SLisandro Dalcin 507818efac9SLisandro Dalcin Level: beginner 508818efac9SLisandro Dalcin 509818efac9SLisandro Dalcin References: 510818efac9SLisandro Dalcin J. Chung, G.M.Hubert. "A Time Integration Algorithm for Structural 511818efac9SLisandro Dalcin Dynamics with Improved Numerical Dissipation: The Generalized-alpha 512818efac9SLisandro Dalcin Method" ASME Journal of Applied Mechanics, 60, 371:375, 1993. 513818efac9SLisandro Dalcin 514818efac9SLisandro Dalcin .seealso: TS, TSCreate(), TSSetType(), TSAlpha2SetRadius(), TSAlpha2SetParams() 515818efac9SLisandro Dalcin M*/ 516818efac9SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSCreate_Alpha2(TS ts) 517818efac9SLisandro Dalcin { 518818efac9SLisandro Dalcin TS_Alpha *th; 519818efac9SLisandro Dalcin PetscErrorCode ierr; 520818efac9SLisandro Dalcin 521818efac9SLisandro Dalcin PetscFunctionBegin; 522818efac9SLisandro Dalcin ts->ops->reset = TSReset_Alpha; 523818efac9SLisandro Dalcin ts->ops->destroy = TSDestroy_Alpha; 524818efac9SLisandro Dalcin ts->ops->view = TSView_Alpha; 525818efac9SLisandro Dalcin ts->ops->setup = TSSetUp_Alpha; 526818efac9SLisandro Dalcin ts->ops->setfromoptions = TSSetFromOptions_Alpha; 527818efac9SLisandro Dalcin ts->ops->step = TSStep_Alpha; 528818efac9SLisandro Dalcin ts->ops->evaluatewlte = TSEvaluateWLTE_Alpha; 529818efac9SLisandro Dalcin ts->ops->rollback = TSRollBack_Alpha; 530818efac9SLisandro Dalcin /*ts->ops->interpolate = TSInterpolate_Alpha;*/ 531818efac9SLisandro Dalcin ts->ops->snesfunction = SNESTSFormFunction_Alpha; 532818efac9SLisandro Dalcin ts->ops->snesjacobian = SNESTSFormJacobian_Alpha; 5332ffb9264SLisandro Dalcin ts->default_adapt_type = TSADAPTNONE; 534818efac9SLisandro Dalcin 535*efd4aadfSBarry Smith ts->usessnes = PETSC_TRUE; 536*efd4aadfSBarry Smith 537818efac9SLisandro Dalcin ierr = PetscNewLog(ts,&th);CHKERRQ(ierr); 538818efac9SLisandro Dalcin ts->data = (void*)th; 539818efac9SLisandro Dalcin 540818efac9SLisandro Dalcin th->Alpha_m = 0.5; 541818efac9SLisandro Dalcin th->Alpha_f = 0.5; 542818efac9SLisandro Dalcin th->Gamma = 0.5; 543818efac9SLisandro Dalcin th->Beta = 0.25; 544818efac9SLisandro Dalcin th->order = 2; 545818efac9SLisandro Dalcin 546818efac9SLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)ts,"TSAlpha2SetRadius_C",TSAlpha2SetRadius_Alpha);CHKERRQ(ierr); 547818efac9SLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)ts,"TSAlpha2SetParams_C",TSAlpha2SetParams_Alpha);CHKERRQ(ierr); 548818efac9SLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)ts,"TSAlpha2GetParams_C",TSAlpha2GetParams_Alpha);CHKERRQ(ierr); 549818efac9SLisandro Dalcin PetscFunctionReturn(0); 550818efac9SLisandro Dalcin } 551818efac9SLisandro Dalcin 552818efac9SLisandro Dalcin /*@ 553818efac9SLisandro Dalcin TSAlpha2SetRadius - sets the desired spectral radius of the method 554818efac9SLisandro Dalcin (i.e. high-frequency numerical damping) 555818efac9SLisandro Dalcin 556818efac9SLisandro Dalcin Logically Collective on TS 557818efac9SLisandro Dalcin 558818efac9SLisandro Dalcin The algorithmic parameters \alpha_m and \alpha_f of the 559818efac9SLisandro Dalcin generalized-\alpha method can be computed in terms of a specified 560818efac9SLisandro Dalcin spectral radius \rho in [0,1] for infinite time step in order to 561818efac9SLisandro Dalcin control high-frequency numerical damping: 562818efac9SLisandro Dalcin \alpha_m = (2-\rho)/(1+\rho) 563818efac9SLisandro Dalcin \alpha_f = 1/(1+\rho) 564818efac9SLisandro Dalcin 565818efac9SLisandro Dalcin Input Parameter: 566818efac9SLisandro Dalcin + ts - timestepping context 567818efac9SLisandro Dalcin - radius - the desired spectral radius 568818efac9SLisandro Dalcin 569818efac9SLisandro Dalcin Options Database: 570818efac9SLisandro Dalcin . -ts_alpha_radius <radius> 571818efac9SLisandro Dalcin 572818efac9SLisandro Dalcin Level: intermediate 573818efac9SLisandro Dalcin 574818efac9SLisandro Dalcin .seealso: TSAlpha2SetParams(), TSAlpha2GetParams() 575818efac9SLisandro Dalcin @*/ 576818efac9SLisandro Dalcin PetscErrorCode TSAlpha2SetRadius(TS ts,PetscReal radius) 577818efac9SLisandro Dalcin { 578818efac9SLisandro Dalcin PetscErrorCode ierr; 579818efac9SLisandro Dalcin 580818efac9SLisandro Dalcin PetscFunctionBegin; 581818efac9SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 582818efac9SLisandro Dalcin PetscValidLogicalCollectiveReal(ts,radius,2); 583818efac9SLisandro Dalcin if (radius < 0 || radius > 1) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Radius %g not in range [0,1]",(double)radius); 584818efac9SLisandro Dalcin ierr = PetscTryMethod(ts,"TSAlpha2SetRadius_C",(TS,PetscReal),(ts,radius));CHKERRQ(ierr); 585818efac9SLisandro Dalcin PetscFunctionReturn(0); 586818efac9SLisandro Dalcin } 587818efac9SLisandro Dalcin 588818efac9SLisandro Dalcin /*@ 589818efac9SLisandro Dalcin TSAlpha2SetParams - sets the algorithmic parameters for TSALPHA2 590818efac9SLisandro Dalcin 591818efac9SLisandro Dalcin Logically Collective on TS 592818efac9SLisandro Dalcin 593818efac9SLisandro Dalcin Second-order accuracy can be obtained so long as: 594818efac9SLisandro Dalcin \gamma = 1/2 + alpha_m - alpha_f 595818efac9SLisandro Dalcin \beta = 1/4 (1 + alpha_m - alpha_f)^2 596818efac9SLisandro Dalcin 597818efac9SLisandro Dalcin Unconditional stability requires: 598818efac9SLisandro Dalcin \alpha_m >= \alpha_f >= 1/2 599818efac9SLisandro Dalcin 600818efac9SLisandro Dalcin 601818efac9SLisandro Dalcin Input Parameter: 602818efac9SLisandro Dalcin + ts - timestepping context 603818efac9SLisandro Dalcin . \alpha_m - algorithmic paramenter 604818efac9SLisandro Dalcin . \alpha_f - algorithmic paramenter 605818efac9SLisandro Dalcin . \gamma - algorithmic paramenter 606818efac9SLisandro Dalcin - \beta - algorithmic paramenter 607818efac9SLisandro Dalcin 608818efac9SLisandro Dalcin Options Database: 609818efac9SLisandro Dalcin + -ts_alpha_alpha_m <alpha_m> 610818efac9SLisandro Dalcin . -ts_alpha_alpha_f <alpha_f> 611818efac9SLisandro Dalcin . -ts_alpha_gamma <gamma> 612818efac9SLisandro Dalcin - -ts_alpha_beta <beta> 613818efac9SLisandro Dalcin 614818efac9SLisandro Dalcin Note: 615818efac9SLisandro Dalcin Use of this function is normally only required to hack TSALPHA2 to 616818efac9SLisandro Dalcin use a modified integration scheme. Users should call 617818efac9SLisandro Dalcin TSAlpha2SetRadius() to set the desired spectral radius of the methods 618818efac9SLisandro Dalcin (i.e. high-frequency damping) in order so select optimal values for 619818efac9SLisandro Dalcin these parameters. 620818efac9SLisandro Dalcin 621818efac9SLisandro Dalcin Level: advanced 622818efac9SLisandro Dalcin 623818efac9SLisandro Dalcin .seealso: TSAlpha2SetRadius(), TSAlpha2GetParams() 624818efac9SLisandro Dalcin @*/ 625818efac9SLisandro Dalcin PetscErrorCode TSAlpha2SetParams(TS ts,PetscReal alpha_m,PetscReal alpha_f,PetscReal gamma,PetscReal beta) 626818efac9SLisandro Dalcin { 627818efac9SLisandro Dalcin PetscErrorCode ierr; 628818efac9SLisandro Dalcin 629818efac9SLisandro Dalcin PetscFunctionBegin; 630818efac9SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 631818efac9SLisandro Dalcin PetscValidLogicalCollectiveReal(ts,alpha_m,2); 632818efac9SLisandro Dalcin PetscValidLogicalCollectiveReal(ts,alpha_f,3); 633818efac9SLisandro Dalcin PetscValidLogicalCollectiveReal(ts,gamma,4); 634818efac9SLisandro Dalcin PetscValidLogicalCollectiveReal(ts,beta,5); 635818efac9SLisandro Dalcin ierr = PetscTryMethod(ts,"TSAlpha2SetParams_C",(TS,PetscReal,PetscReal,PetscReal,PetscReal),(ts,alpha_m,alpha_f,gamma,beta));CHKERRQ(ierr); 636818efac9SLisandro Dalcin PetscFunctionReturn(0); 637818efac9SLisandro Dalcin } 638818efac9SLisandro Dalcin 639818efac9SLisandro Dalcin /*@ 640818efac9SLisandro Dalcin TSAlpha2GetParams - gets the algorithmic parameters for TSALPHA2 641818efac9SLisandro Dalcin 642818efac9SLisandro Dalcin Not Collective 643818efac9SLisandro Dalcin 644818efac9SLisandro Dalcin Input Parameter: 645818efac9SLisandro Dalcin . ts - timestepping context 646818efac9SLisandro Dalcin 647818efac9SLisandro Dalcin Output Parameters: 648818efac9SLisandro Dalcin + \alpha_m - algorithmic parameter 649818efac9SLisandro Dalcin . \alpha_f - algorithmic parameter 650818efac9SLisandro Dalcin . \gamma - algorithmic parameter 651818efac9SLisandro Dalcin - \beta - algorithmic parameter 652818efac9SLisandro Dalcin 653818efac9SLisandro Dalcin Note: 654818efac9SLisandro Dalcin Use of this function is normally only required to hack TSALPHA2 to 655818efac9SLisandro Dalcin use a modified integration scheme. Users should call 656818efac9SLisandro Dalcin TSAlpha2SetRadius() to set the high-frequency damping (i.e. spectral 657818efac9SLisandro Dalcin radius of the method) in order so select optimal values for these 658818efac9SLisandro Dalcin parameters. 659818efac9SLisandro Dalcin 660818efac9SLisandro Dalcin Level: advanced 661818efac9SLisandro Dalcin 662818efac9SLisandro Dalcin .seealso: TSAlpha2SetRadius(), TSAlpha2SetParams() 663818efac9SLisandro Dalcin @*/ 664818efac9SLisandro Dalcin PetscErrorCode TSAlpha2GetParams(TS ts,PetscReal *alpha_m,PetscReal *alpha_f,PetscReal *gamma,PetscReal *beta) 665818efac9SLisandro Dalcin { 666818efac9SLisandro Dalcin PetscErrorCode ierr; 667818efac9SLisandro Dalcin 668818efac9SLisandro Dalcin PetscFunctionBegin; 669818efac9SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 670818efac9SLisandro Dalcin if (alpha_m) PetscValidRealPointer(alpha_m,2); 671818efac9SLisandro Dalcin if (alpha_f) PetscValidRealPointer(alpha_f,3); 672818efac9SLisandro Dalcin if (gamma) PetscValidRealPointer(gamma,4); 673818efac9SLisandro Dalcin if (beta) PetscValidRealPointer(beta,5); 674818efac9SLisandro Dalcin ierr = PetscUseMethod(ts,"TSAlpha2GetParams_C",(TS,PetscReal*,PetscReal*,PetscReal*,PetscReal*),(ts,alpha_m,alpha_f,gamma,beta));CHKERRQ(ierr); 675818efac9SLisandro Dalcin PetscFunctionReturn(0); 676818efac9SLisandro Dalcin } 677