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 PetscBool adapt; 33b07a2398SLisandro Dalcin Vec vec_sol_prev; 341566a47fSLisandro Dalcin Vec vec_lte_work; 35b07a2398SLisandro Dalcin 36b07a2398SLisandro Dalcin TSStepStatus status; 37b07a2398SLisandro Dalcin } TS_Alpha; 38b07a2398SLisandro Dalcin 39b07a2398SLisandro Dalcin #undef __FUNCT__ 40b07a2398SLisandro Dalcin #define __FUNCT__ "TSAlpha_StageTime" 41b07a2398SLisandro Dalcin static PetscErrorCode TSAlpha_StageTime(TS ts) 42b07a2398SLisandro Dalcin { 43b07a2398SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 44b07a2398SLisandro Dalcin PetscReal t = ts->ptime; 45b07a2398SLisandro Dalcin PetscReal dt = ts->time_step; 46b07a2398SLisandro Dalcin PetscReal Alpha_m = th->Alpha_m; 47b07a2398SLisandro Dalcin PetscReal Alpha_f = th->Alpha_f; 48b07a2398SLisandro Dalcin PetscReal Gamma = th->Gamma; 49b07a2398SLisandro Dalcin 50b07a2398SLisandro Dalcin PetscFunctionBegin; 51b07a2398SLisandro Dalcin th->stage_time = t + Alpha_f*dt; 52b07a2398SLisandro Dalcin th->shift_V = Alpha_m/(Alpha_f*Gamma*dt); 53b07a2398SLisandro Dalcin th->scale_F = 1/Alpha_f; 54b07a2398SLisandro Dalcin PetscFunctionReturn(0); 55b07a2398SLisandro Dalcin } 56b07a2398SLisandro Dalcin 57b07a2398SLisandro Dalcin #undef __FUNCT__ 58b07a2398SLisandro Dalcin #define __FUNCT__ "TSAlpha_StageVecs" 59b07a2398SLisandro Dalcin static PetscErrorCode TSAlpha_StageVecs(TS ts,Vec X) 60b07a2398SLisandro Dalcin { 61b07a2398SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 62b07a2398SLisandro Dalcin Vec X1 = X, V1 = th->V1; 63b07a2398SLisandro Dalcin Vec Xa = th->Xa, Va = th->Va; 64b07a2398SLisandro Dalcin Vec X0 = th->X0, V0 = th->V0; 65b07a2398SLisandro Dalcin PetscReal dt = ts->time_step; 66b07a2398SLisandro Dalcin PetscReal Alpha_m = th->Alpha_m; 67b07a2398SLisandro Dalcin PetscReal Alpha_f = th->Alpha_f; 68b07a2398SLisandro Dalcin PetscReal Gamma = th->Gamma; 69b07a2398SLisandro Dalcin PetscErrorCode ierr; 70b07a2398SLisandro Dalcin 71b07a2398SLisandro Dalcin PetscFunctionBegin; 72b07a2398SLisandro Dalcin /* V1 = 1/(Gamma*dT)*(X1-X0) + (1-1/Gamma)*V0 */ 73b07a2398SLisandro Dalcin ierr = VecWAXPY(V1,-1.0,X0,X1);CHKERRQ(ierr); 74b07a2398SLisandro Dalcin ierr = VecAXPBY(V1,1-1/Gamma,1/(Gamma*dt),V0);CHKERRQ(ierr); 75b07a2398SLisandro Dalcin /* Xa = X0 + Alpha_f*(X1-X0) */ 76b07a2398SLisandro Dalcin ierr = VecWAXPY(Xa,-1.0,X0,X1);CHKERRQ(ierr); 77b07a2398SLisandro Dalcin ierr = VecAYPX(Xa,Alpha_f,X0);CHKERRQ(ierr); 78b07a2398SLisandro Dalcin /* Va = V0 + Alpha_m*(V1-V0) */ 79b07a2398SLisandro Dalcin ierr = VecWAXPY(Va,-1.0,V0,V1);CHKERRQ(ierr); 80b07a2398SLisandro Dalcin ierr = VecAYPX(Va,Alpha_m,V0);CHKERRQ(ierr); 81b07a2398SLisandro Dalcin PetscFunctionReturn(0); 82b07a2398SLisandro Dalcin } 83b07a2398SLisandro Dalcin 84b07a2398SLisandro Dalcin #undef __FUNCT__ 85b07a2398SLisandro Dalcin #define __FUNCT__ "TS_SNESSolve" 86b07a2398SLisandro Dalcin static PetscErrorCode TS_SNESSolve(TS ts,Vec b,Vec x) 87b07a2398SLisandro Dalcin { 88b07a2398SLisandro Dalcin PetscInt nits,lits; 89b07a2398SLisandro Dalcin PetscErrorCode ierr; 90b07a2398SLisandro Dalcin 91b07a2398SLisandro Dalcin PetscFunctionBegin; 92b07a2398SLisandro Dalcin ierr = SNESSolve(ts->snes,b,x);CHKERRQ(ierr); 93b07a2398SLisandro Dalcin ierr = SNESGetIterationNumber(ts->snes,&nits);CHKERRQ(ierr); 94b07a2398SLisandro Dalcin ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr); 95b07a2398SLisandro Dalcin ts->snes_its += nits; ts->ksp_its += lits; 96b07a2398SLisandro Dalcin PetscFunctionReturn(0); 97b07a2398SLisandro Dalcin } 98b07a2398SLisandro Dalcin 99b07a2398SLisandro Dalcin /* 100b07a2398SLisandro Dalcin Compute a consistent initial state for the generalized-alpha method. 101b07a2398SLisandro Dalcin - Solve two successive backward Euler steps with halved time step. 102b07a2398SLisandro Dalcin - Compute the initial time derivative using backward differences. 103b07a2398SLisandro Dalcin - If using adaptivity, estimate the LTE of the initial step. 104b07a2398SLisandro Dalcin */ 105b07a2398SLisandro Dalcin #undef __FUNCT__ 106b07a2398SLisandro Dalcin #define __FUNCT__ "TSAlpha_ResetStep" 107b07a2398SLisandro Dalcin static PetscErrorCode TSAlpha_ResetStep(TS ts,PetscBool *initok) 108b07a2398SLisandro Dalcin { 109b07a2398SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 110b07a2398SLisandro Dalcin PetscReal time_step; 111b07a2398SLisandro Dalcin PetscReal alpha_m,alpha_f,gamma; 112b07a2398SLisandro Dalcin Vec X0 = ts->vec_sol, X1, X2 = th->X1; 113b07a2398SLisandro Dalcin PetscBool stageok; 114b07a2398SLisandro Dalcin PetscErrorCode ierr; 115b07a2398SLisandro Dalcin 116b07a2398SLisandro Dalcin PetscFunctionBegin; 117b07a2398SLisandro Dalcin ierr = VecDuplicate(X0,&X1);CHKERRQ(ierr); 118b07a2398SLisandro Dalcin 119b07a2398SLisandro Dalcin /* Setup backward Euler with halved time step */ 120b07a2398SLisandro Dalcin ierr = TSAlphaGetParams(ts,&alpha_m,&alpha_f,&gamma);CHKERRQ(ierr); 121b07a2398SLisandro Dalcin ierr = TSAlphaSetParams(ts,1,1,1);CHKERRQ(ierr); 122b07a2398SLisandro Dalcin ierr = TSGetTimeStep(ts,&time_step);CHKERRQ(ierr); 123b07a2398SLisandro Dalcin ts->time_step = time_step/2; 124b07a2398SLisandro Dalcin ierr = TSAlpha_StageTime(ts);CHKERRQ(ierr); 125b07a2398SLisandro Dalcin th->stage_time = ts->ptime; 126b07a2398SLisandro Dalcin ierr = VecZeroEntries(th->V0);CHKERRQ(ierr); 127b07a2398SLisandro Dalcin 128b07a2398SLisandro Dalcin /* First BE step, (t0,X0) -> (t1,X1) */ 129b07a2398SLisandro Dalcin th->stage_time += ts->time_step; 130b07a2398SLisandro Dalcin ierr = VecCopy(X0,th->X0);CHKERRQ(ierr); 131b07a2398SLisandro Dalcin ierr = TSPreStage(ts,th->stage_time);CHKERRQ(ierr); 132b07a2398SLisandro Dalcin ierr = VecCopy(th->X0,X1);CHKERRQ(ierr); 133b07a2398SLisandro Dalcin ierr = TS_SNESSolve(ts,NULL,X1);CHKERRQ(ierr); 134b07a2398SLisandro Dalcin ierr = TSPostStage(ts,th->stage_time,0,&X1);CHKERRQ(ierr); 135b07a2398SLisandro Dalcin ierr = TSAdaptCheckStage(ts->adapt,ts,th->stage_time,X1,&stageok);CHKERRQ(ierr); 136b07a2398SLisandro Dalcin if (!stageok) goto finally; 137b07a2398SLisandro Dalcin 138b07a2398SLisandro Dalcin /* Second BE step, (t1,X1) -> (t2,X2) */ 139b07a2398SLisandro Dalcin th->stage_time += ts->time_step; 140b07a2398SLisandro Dalcin ierr = VecCopy(X1,th->X0);CHKERRQ(ierr); 141b07a2398SLisandro Dalcin ierr = TSPreStage(ts,th->stage_time);CHKERRQ(ierr); 142b07a2398SLisandro Dalcin ierr = VecCopy(th->X0,X2);CHKERRQ(ierr); 143b07a2398SLisandro Dalcin ierr = TS_SNESSolve(ts,NULL,X2);CHKERRQ(ierr); 144b07a2398SLisandro Dalcin ierr = TSPostStage(ts,th->stage_time,0,&X2);CHKERRQ(ierr); 145b07a2398SLisandro Dalcin ierr = TSAdaptCheckStage(ts->adapt,ts,th->stage_time,X2,&stageok);CHKERRQ(ierr); 146b07a2398SLisandro Dalcin if (!stageok) goto finally; 147b07a2398SLisandro Dalcin 148b07a2398SLisandro Dalcin /* Compute V0 ~ dX/dt at t0 with backward differences */ 149b07a2398SLisandro Dalcin ierr = VecZeroEntries(th->V0);CHKERRQ(ierr); 150b07a2398SLisandro Dalcin ierr = VecAXPY(th->V0,-3/ts->time_step,X0);CHKERRQ(ierr); 151b07a2398SLisandro Dalcin ierr = VecAXPY(th->V0,+4/ts->time_step,X1);CHKERRQ(ierr); 152b07a2398SLisandro Dalcin ierr = VecAXPY(th->V0,-1/ts->time_step,X2);CHKERRQ(ierr); 153b07a2398SLisandro Dalcin 154b07a2398SLisandro Dalcin /* Rough, lower-order estimate LTE of the initial step */ 155b07a2398SLisandro Dalcin if (th->adapt) { 156b07a2398SLisandro Dalcin ierr = VecZeroEntries(th->vec_sol_prev);CHKERRQ(ierr); 157b07a2398SLisandro Dalcin ierr = VecAXPY(th->vec_sol_prev,+2,X2);CHKERRQ(ierr); 158b07a2398SLisandro Dalcin ierr = VecAXPY(th->vec_sol_prev,-4,X1);CHKERRQ(ierr); 159b07a2398SLisandro Dalcin ierr = VecAXPY(th->vec_sol_prev,+2,X0);CHKERRQ(ierr); 160b07a2398SLisandro Dalcin } 161b07a2398SLisandro Dalcin 162b07a2398SLisandro Dalcin finally: 163b07a2398SLisandro Dalcin /* Revert TSAlpha to the initial state (t0,X0) */ 164b07a2398SLisandro Dalcin if (initok) *initok = stageok; 165b07a2398SLisandro Dalcin ierr = TSSetTimeStep(ts,time_step);CHKERRQ(ierr); 166b07a2398SLisandro Dalcin ierr = TSAlphaSetParams(ts,alpha_m,alpha_f,gamma);CHKERRQ(ierr); 167b07a2398SLisandro Dalcin ierr = VecCopy(ts->vec_sol,th->X0);CHKERRQ(ierr); 168b07a2398SLisandro Dalcin 169b07a2398SLisandro Dalcin ierr = VecDestroy(&X1);CHKERRQ(ierr); 170b07a2398SLisandro Dalcin PetscFunctionReturn(0); 171b07a2398SLisandro Dalcin } 172b07a2398SLisandro Dalcin 173b07a2398SLisandro Dalcin #define TSEvent_Status(ts) (ts->event ? ts->event->status : TSEVENT_NONE) 174b07a2398SLisandro Dalcin 175b07a2398SLisandro Dalcin #undef __FUNCT__ 176b07a2398SLisandro Dalcin #define __FUNCT__ "TSStep_Alpha" 177b07a2398SLisandro Dalcin static PetscErrorCode TSStep_Alpha(TS ts) 178b07a2398SLisandro Dalcin { 179b07a2398SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 180b07a2398SLisandro Dalcin PetscInt rejections = 0; 181b07a2398SLisandro Dalcin PetscBool stageok,accept = PETSC_TRUE; 182b07a2398SLisandro Dalcin PetscReal next_time_step = ts->time_step; 183b07a2398SLisandro Dalcin PetscErrorCode ierr; 184b07a2398SLisandro Dalcin 185b07a2398SLisandro Dalcin PetscFunctionBegin; 186b07a2398SLisandro Dalcin ierr = PetscCitationsRegister(citation,&cited);CHKERRQ(ierr); 187b07a2398SLisandro Dalcin 188b07a2398SLisandro Dalcin if (!ts->steprollback) { 189b07a2398SLisandro Dalcin if (th->adapt) { ierr = VecCopy(th->X0,th->vec_sol_prev);CHKERRQ(ierr); } 190b07a2398SLisandro Dalcin ierr = VecCopy(ts->vec_sol,th->X0);CHKERRQ(ierr); 191b07a2398SLisandro Dalcin ierr = VecCopy(th->V1,th->V0);CHKERRQ(ierr); 192b07a2398SLisandro Dalcin } 193b07a2398SLisandro Dalcin 1941566a47fSLisandro Dalcin th->status = TS_STEP_INCOMPLETE; 195b07a2398SLisandro Dalcin while (!ts->reason && th->status != TS_STEP_COMPLETE) { 196b07a2398SLisandro Dalcin 197b07a2398SLisandro Dalcin if (!ts->steps || TSEvent_Status(ts) == TSEVENT_RESET_NEXTSTEP) { 198b07a2398SLisandro Dalcin ierr = TSAlpha_ResetStep(ts,&stageok);CHKERRQ(ierr); 199b07a2398SLisandro Dalcin if (!stageok) {accept = PETSC_FALSE; goto reject_step;} 200b07a2398SLisandro Dalcin } 201b07a2398SLisandro Dalcin 202b07a2398SLisandro Dalcin ierr = TSAlpha_StageTime(ts);CHKERRQ(ierr); 203b07a2398SLisandro Dalcin ierr = VecCopy(th->X0,th->X1);CHKERRQ(ierr); 204b07a2398SLisandro Dalcin ierr = TSPreStage(ts,th->stage_time);CHKERRQ(ierr); 205b07a2398SLisandro Dalcin ierr = TS_SNESSolve(ts,NULL,th->X1);CHKERRQ(ierr); 206*be5899b3SLisandro Dalcin ierr = TSPostStage(ts,th->stage_time,0,&th->Xa);CHKERRQ(ierr); 207*be5899b3SLisandro Dalcin ierr = TSAdaptCheckStage(ts->adapt,ts,th->stage_time,th->Xa,&stageok);CHKERRQ(ierr); 208b07a2398SLisandro Dalcin if (!stageok) {accept = PETSC_FALSE; goto reject_step;} 209b07a2398SLisandro Dalcin 2101566a47fSLisandro Dalcin th->status = TS_STEP_PENDING; 211b07a2398SLisandro Dalcin ierr = VecCopy(th->X1,ts->vec_sol);CHKERRQ(ierr); 2121566a47fSLisandro Dalcin ierr = TSAdaptChoose(ts->adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr); 2131566a47fSLisandro Dalcin th->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 214b07a2398SLisandro Dalcin if (!accept) { 215b07a2398SLisandro Dalcin ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr); 216*be5899b3SLisandro Dalcin ts->time_step = next_time_step; 217*be5899b3SLisandro Dalcin goto reject_step; 218b07a2398SLisandro Dalcin } 219b07a2398SLisandro Dalcin 220b07a2398SLisandro Dalcin ts->ptime += ts->time_step; 221b07a2398SLisandro Dalcin ts->time_step = next_time_step; 222b07a2398SLisandro Dalcin break; 223b07a2398SLisandro Dalcin 224b07a2398SLisandro Dalcin reject_step: 225b07a2398SLisandro Dalcin ts->reject++; 226b07a2398SLisandro Dalcin if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 227b07a2398SLisandro Dalcin ts->reason = TS_DIVERGED_STEP_REJECTED; 228b07a2398SLisandro Dalcin ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr); 229b07a2398SLisandro Dalcin } 2301566a47fSLisandro Dalcin 231b07a2398SLisandro Dalcin } 232b07a2398SLisandro Dalcin PetscFunctionReturn(0); 233b07a2398SLisandro Dalcin } 234b07a2398SLisandro Dalcin 235b07a2398SLisandro Dalcin #undef __FUNCT__ 2369808bdc1SLisandro Dalcin #define __FUNCT__ "TSEvaluateWLTE_Alpha" 2379808bdc1SLisandro Dalcin static PetscErrorCode TSEvaluateWLTE_Alpha(TS ts,NormType wnormtype,PetscInt *order,PetscReal *wlte) 238b07a2398SLisandro Dalcin { 239b07a2398SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 2409808bdc1SLisandro Dalcin Vec X = th->X1; /* X = solution */ 2411566a47fSLisandro Dalcin Vec Y = th->vec_lte_work; /* Y = X + LTE */ 242b07a2398SLisandro Dalcin PetscErrorCode ierr; 243b07a2398SLisandro Dalcin 244b07a2398SLisandro Dalcin PetscFunctionBegin; 245b07a2398SLisandro Dalcin if (!ts->steps || TSEvent_Status(ts) == TSEVENT_RESET_NEXTSTEP) { 246057721c2SLisandro Dalcin /* th->vec_sol_prev is set to the LTE in TSAlpha_ResetStep() */ 2479808bdc1SLisandro Dalcin ierr = VecWAXPY(Y,1.0,th->vec_sol_prev,X);CHKERRQ(ierr); 248b07a2398SLisandro Dalcin } else { 249b07a2398SLisandro Dalcin /* Compute LTE using backward differences with non-constant time step */ 250*be5899b3SLisandro Dalcin PetscReal h = ts->time_step, h_prev = ts->ptime - ts->ptime_prev; 251*be5899b3SLisandro Dalcin PetscReal a = 1 + h_prev/h; 252b07a2398SLisandro Dalcin PetscScalar scal[3]; Vec vecs[3]; 253b07a2398SLisandro Dalcin scal[0] = +1/a; scal[1] = -1/(a-1); scal[2] = +1/(a*(a-1)); 254b07a2398SLisandro Dalcin vecs[0] = th->X1; vecs[1] = th->X0; vecs[2] = th->vec_sol_prev; 2559808bdc1SLisandro Dalcin ierr = VecCopy(X,Y);CHKERRQ(ierr); 2569808bdc1SLisandro Dalcin ierr = VecMAXPY(Y,3,scal,vecs);CHKERRQ(ierr); 257b07a2398SLisandro Dalcin } 2589808bdc1SLisandro Dalcin ierr = TSErrorWeightedNorm(ts,X,Y,wnormtype,wlte);CHKERRQ(ierr); 2599808bdc1SLisandro Dalcin if (order) *order = 2; 260b07a2398SLisandro Dalcin PetscFunctionReturn(0); 261b07a2398SLisandro Dalcin } 262b07a2398SLisandro Dalcin 263b07a2398SLisandro Dalcin #undef __FUNCT__ 264b07a2398SLisandro Dalcin #define __FUNCT__ "TSRollBack_Alpha" 265b07a2398SLisandro Dalcin static PetscErrorCode TSRollBack_Alpha(TS ts) 266b07a2398SLisandro Dalcin { 267b07a2398SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 268b07a2398SLisandro Dalcin PetscErrorCode ierr; 269b07a2398SLisandro Dalcin 270b07a2398SLisandro Dalcin PetscFunctionBegin; 271b07a2398SLisandro Dalcin ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr); 272b07a2398SLisandro Dalcin PetscFunctionReturn(0); 273b07a2398SLisandro Dalcin } 274b07a2398SLisandro Dalcin 275b07a2398SLisandro Dalcin #undef __FUNCT__ 276b07a2398SLisandro Dalcin #define __FUNCT__ "TSInterpolate_Alpha" 277b07a2398SLisandro Dalcin static PetscErrorCode TSInterpolate_Alpha(TS ts,PetscReal t,Vec X) 278b07a2398SLisandro Dalcin { 279b07a2398SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 280b07a2398SLisandro Dalcin PetscReal dt = t - ts->ptime; 281b07a2398SLisandro Dalcin PetscErrorCode ierr; 282b07a2398SLisandro Dalcin 283b07a2398SLisandro Dalcin PetscFunctionBegin; 284b07a2398SLisandro Dalcin ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 285b07a2398SLisandro Dalcin ierr = VecAXPY(X,th->Gamma*dt,th->V1);CHKERRQ(ierr); 286b07a2398SLisandro Dalcin ierr = VecAXPY(X,(1-th->Gamma)*dt,th->V0);CHKERRQ(ierr); 287b07a2398SLisandro Dalcin PetscFunctionReturn(0); 288b07a2398SLisandro Dalcin } 289b07a2398SLisandro Dalcin 290b07a2398SLisandro Dalcin #undef __FUNCT__ 291b07a2398SLisandro Dalcin #define __FUNCT__ "SNESTSFormFunction_Alpha" 292b07a2398SLisandro Dalcin static PetscErrorCode SNESTSFormFunction_Alpha(PETSC_UNUSED SNES snes,Vec X,Vec F,TS ts) 293b07a2398SLisandro Dalcin { 294b07a2398SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 295b07a2398SLisandro Dalcin PetscReal ta = th->stage_time; 296b07a2398SLisandro Dalcin Vec Xa = th->Xa, Va = th->Va; 297b07a2398SLisandro Dalcin PetscErrorCode ierr; 298b07a2398SLisandro Dalcin 299b07a2398SLisandro Dalcin PetscFunctionBegin; 300b07a2398SLisandro Dalcin ierr = TSAlpha_StageVecs(ts,X);CHKERRQ(ierr); 301b07a2398SLisandro Dalcin /* F = Function(ta,Xa,Va) */ 302b07a2398SLisandro Dalcin ierr = TSComputeIFunction(ts,ta,Xa,Va,F,PETSC_FALSE);CHKERRQ(ierr); 303b07a2398SLisandro Dalcin ierr = VecScale(F,th->scale_F);CHKERRQ(ierr); 304b07a2398SLisandro Dalcin PetscFunctionReturn(0); 305b07a2398SLisandro Dalcin } 306b07a2398SLisandro Dalcin 307b07a2398SLisandro Dalcin #undef __FUNCT__ 308b07a2398SLisandro Dalcin #define __FUNCT__ "SNESTSFormJacobian_Alpha" 309b07a2398SLisandro Dalcin static PetscErrorCode SNESTSFormJacobian_Alpha(PETSC_UNUSED SNES snes,PETSC_UNUSED Vec X,Mat J,Mat P,TS ts) 310b07a2398SLisandro Dalcin { 311b07a2398SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 312b07a2398SLisandro Dalcin PetscReal ta = th->stage_time; 313b07a2398SLisandro Dalcin Vec Xa = th->Xa, Va = th->Va; 314b07a2398SLisandro Dalcin PetscReal dVdX = th->shift_V; 315b07a2398SLisandro Dalcin PetscErrorCode ierr; 316b07a2398SLisandro Dalcin 317b07a2398SLisandro Dalcin PetscFunctionBegin; 318b07a2398SLisandro Dalcin /* J,P = Jacobian(ta,Xa,Va) */ 319b07a2398SLisandro Dalcin ierr = TSComputeIJacobian(ts,ta,Xa,Va,dVdX,J,P,PETSC_FALSE);CHKERRQ(ierr); 320b07a2398SLisandro Dalcin PetscFunctionReturn(0); 321b07a2398SLisandro Dalcin } 322b07a2398SLisandro Dalcin 323b07a2398SLisandro Dalcin #undef __FUNCT__ 324b07a2398SLisandro Dalcin #define __FUNCT__ "TSReset_Alpha" 325b07a2398SLisandro Dalcin static PetscErrorCode TSReset_Alpha(TS ts) 326b07a2398SLisandro Dalcin { 327b07a2398SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 328b07a2398SLisandro Dalcin PetscErrorCode ierr; 329b07a2398SLisandro Dalcin 330b07a2398SLisandro Dalcin PetscFunctionBegin; 331b07a2398SLisandro Dalcin ierr = VecDestroy(&th->X0);CHKERRQ(ierr); 332b07a2398SLisandro Dalcin ierr = VecDestroy(&th->Xa);CHKERRQ(ierr); 333b07a2398SLisandro Dalcin ierr = VecDestroy(&th->X1);CHKERRQ(ierr); 334b07a2398SLisandro Dalcin ierr = VecDestroy(&th->V0);CHKERRQ(ierr); 335b07a2398SLisandro Dalcin ierr = VecDestroy(&th->Va);CHKERRQ(ierr); 336b07a2398SLisandro Dalcin ierr = VecDestroy(&th->V1);CHKERRQ(ierr); 337b07a2398SLisandro Dalcin ierr = VecDestroy(&th->vec_sol_prev);CHKERRQ(ierr); 3381566a47fSLisandro Dalcin ierr = VecDestroy(&th->vec_lte_work);CHKERRQ(ierr); 339b07a2398SLisandro Dalcin PetscFunctionReturn(0); 340b07a2398SLisandro Dalcin } 341b07a2398SLisandro Dalcin 342b07a2398SLisandro Dalcin #undef __FUNCT__ 343b07a2398SLisandro Dalcin #define __FUNCT__ "TSDestroy_Alpha" 344b07a2398SLisandro Dalcin static PetscErrorCode TSDestroy_Alpha(TS ts) 345b07a2398SLisandro Dalcin { 346b07a2398SLisandro Dalcin PetscErrorCode ierr; 347b07a2398SLisandro Dalcin 348b07a2398SLisandro Dalcin PetscFunctionBegin; 349b07a2398SLisandro Dalcin ierr = TSReset_Alpha(ts);CHKERRQ(ierr); 350b07a2398SLisandro Dalcin ierr = PetscFree(ts->data);CHKERRQ(ierr); 351b07a2398SLisandro Dalcin 352b07a2398SLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)ts,"TSAlphaUseAdapt_C",NULL);CHKERRQ(ierr); 353b07a2398SLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)ts,"TSAlphaSetRadius_C",NULL);CHKERRQ(ierr); 354b07a2398SLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)ts,"TSAlphaSetParams_C",NULL);CHKERRQ(ierr); 355b07a2398SLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)ts,"TSAlphaGetParams_C",NULL);CHKERRQ(ierr); 356b07a2398SLisandro Dalcin PetscFunctionReturn(0); 357b07a2398SLisandro Dalcin } 358b07a2398SLisandro Dalcin 359b07a2398SLisandro Dalcin #undef __FUNCT__ 360b07a2398SLisandro Dalcin #define __FUNCT__ "TSSetUp_Alpha" 361b07a2398SLisandro Dalcin static PetscErrorCode TSSetUp_Alpha(TS ts) 362b07a2398SLisandro Dalcin { 363b07a2398SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 364b07a2398SLisandro Dalcin PetscErrorCode ierr; 365b07a2398SLisandro Dalcin 366b07a2398SLisandro Dalcin PetscFunctionBegin; 367b07a2398SLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->X0);CHKERRQ(ierr); 368b07a2398SLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->Xa);CHKERRQ(ierr); 369b07a2398SLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->X1);CHKERRQ(ierr); 370b07a2398SLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->V0);CHKERRQ(ierr); 371b07a2398SLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->Va);CHKERRQ(ierr); 372b07a2398SLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->V1);CHKERRQ(ierr); 3731566a47fSLisandro Dalcin 374b07a2398SLisandro Dalcin ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr); 3751566a47fSLisandro Dalcin ierr = TSAdaptCandidatesClear(ts->adapt);CHKERRQ(ierr); 3761566a47fSLisandro Dalcin if (!th->adapt) { 377b07a2398SLisandro Dalcin ierr = TSAdaptSetType(ts->adapt,TSADAPTNONE);CHKERRQ(ierr); 378b07a2398SLisandro Dalcin } else { 379b07a2398SLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->vec_sol_prev);CHKERRQ(ierr); 3801566a47fSLisandro Dalcin ierr = VecDuplicate(ts->vec_sol,&th->vec_lte_work);CHKERRQ(ierr); 381b07a2398SLisandro Dalcin if (ts->exact_final_time == TS_EXACTFINALTIME_UNSPECIFIED) 382b07a2398SLisandro Dalcin ts->exact_final_time = TS_EXACTFINALTIME_MATCHSTEP; 383b07a2398SLisandro Dalcin } 3841566a47fSLisandro Dalcin 385b07a2398SLisandro Dalcin ierr = TSGetSNES(ts,&ts->snes);CHKERRQ(ierr); 386b07a2398SLisandro Dalcin PetscFunctionReturn(0); 387b07a2398SLisandro Dalcin } 388b07a2398SLisandro Dalcin 389b07a2398SLisandro Dalcin #undef __FUNCT__ 390b07a2398SLisandro Dalcin #define __FUNCT__ "TSSetFromOptions_Alpha" 391b07a2398SLisandro Dalcin static PetscErrorCode TSSetFromOptions_Alpha(PetscOptionItems *PetscOptionsObject,TS ts) 392b07a2398SLisandro Dalcin { 393b07a2398SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 394b07a2398SLisandro Dalcin PetscErrorCode ierr; 395b07a2398SLisandro Dalcin 396b07a2398SLisandro Dalcin PetscFunctionBegin; 397b07a2398SLisandro Dalcin ierr = PetscOptionsHead(PetscOptionsObject,"Generalized-Alpha ODE solver options");CHKERRQ(ierr); 398b07a2398SLisandro Dalcin { 399b07a2398SLisandro Dalcin PetscBool flg; 400b07a2398SLisandro Dalcin PetscReal radius = 1; 401b07a2398SLisandro Dalcin PetscBool adapt = th->adapt; 402b07a2398SLisandro Dalcin ierr = PetscOptionsReal("-ts_alpha_radius","Spectral radius (high-frequency dissipation)","TSAlphaSetRadius",radius,&radius,&flg);CHKERRQ(ierr); 403b07a2398SLisandro Dalcin if (flg) {ierr = TSAlphaSetRadius(ts,radius);CHKERRQ(ierr);} 404b07a2398SLisandro Dalcin ierr = PetscOptionsReal("-ts_alpha_alpha_m","Algoritmic parameter alpha_m","TSAlphaSetParams",th->Alpha_m,&th->Alpha_m,NULL);CHKERRQ(ierr); 405b07a2398SLisandro Dalcin ierr = PetscOptionsReal("-ts_alpha_alpha_f","Algoritmic parameter alpha_f","TSAlphaSetParams",th->Alpha_f,&th->Alpha_f,NULL);CHKERRQ(ierr); 406b07a2398SLisandro Dalcin ierr = PetscOptionsReal("-ts_alpha_gamma","Algoritmic parameter gamma","TSAlphaSetParams",th->Gamma,&th->Gamma,NULL);CHKERRQ(ierr); 407b07a2398SLisandro Dalcin ierr = TSAlphaSetParams(ts,th->Alpha_m,th->Alpha_f,th->Gamma);CHKERRQ(ierr); 408b07a2398SLisandro Dalcin ierr = PetscOptionsBool("-ts_alpha_adapt","Use time-step adaptivity with the Alpha method","TSAlpha2UseAdapt",adapt,&adapt,&flg);CHKERRQ(ierr); 409b07a2398SLisandro Dalcin if (flg) {ierr = TSAlphaUseAdapt(ts,adapt);CHKERRQ(ierr);} 410b07a2398SLisandro Dalcin } 411b07a2398SLisandro Dalcin ierr = PetscOptionsTail();CHKERRQ(ierr); 412b07a2398SLisandro Dalcin PetscFunctionReturn(0); 413b07a2398SLisandro Dalcin } 414b07a2398SLisandro Dalcin 415b07a2398SLisandro Dalcin #undef __FUNCT__ 416b07a2398SLisandro Dalcin #define __FUNCT__ "TSView_Alpha" 417b07a2398SLisandro Dalcin static PetscErrorCode TSView_Alpha(TS ts,PetscViewer viewer) 418b07a2398SLisandro Dalcin { 419b07a2398SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 4209c334d8fSLisandro Dalcin PetscBool iascii; 421b07a2398SLisandro Dalcin PetscErrorCode ierr; 422b07a2398SLisandro Dalcin 423b07a2398SLisandro Dalcin PetscFunctionBegin; 4249c334d8fSLisandro Dalcin ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 4259c334d8fSLisandro Dalcin if (iascii) { 4269c334d8fSLisandro 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); 4279c334d8fSLisandro Dalcin } 4289c334d8fSLisandro Dalcin if (ts->adapt) {ierr = TSAdaptView(ts->adapt,viewer);CHKERRQ(ierr);} 429b07a2398SLisandro Dalcin if (ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);} 430b07a2398SLisandro Dalcin PetscFunctionReturn(0); 431b07a2398SLisandro Dalcin } 432b07a2398SLisandro Dalcin 433b07a2398SLisandro Dalcin #undef __FUNCT__ 434b07a2398SLisandro Dalcin #define __FUNCT__ "TSAlphaUseAdapt_Alpha" 435b07a2398SLisandro Dalcin static PetscErrorCode TSAlphaUseAdapt_Alpha(TS ts,PetscBool use) 436b07a2398SLisandro Dalcin { 437b07a2398SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 438b07a2398SLisandro Dalcin 439b07a2398SLisandro Dalcin PetscFunctionBegin; 440b07a2398SLisandro Dalcin if (use == th->adapt) PetscFunctionReturn(0); 441b07a2398SLisandro Dalcin if (ts->setupcalled) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ORDER,"Cannot change adaptivity after TSSetUp()"); 442b07a2398SLisandro Dalcin th->adapt = use; 443b07a2398SLisandro Dalcin PetscFunctionReturn(0); 444b07a2398SLisandro Dalcin } 445b07a2398SLisandro Dalcin 446b07a2398SLisandro Dalcin #undef __FUNCT__ 447b07a2398SLisandro Dalcin #define __FUNCT__ "TSAlphaSetRadius_Alpha" 448b07a2398SLisandro Dalcin static PetscErrorCode TSAlphaSetRadius_Alpha(TS ts,PetscReal radius) 449b07a2398SLisandro Dalcin { 450b07a2398SLisandro Dalcin PetscReal alpha_m,alpha_f,gamma; 451b07a2398SLisandro Dalcin PetscErrorCode ierr; 452b07a2398SLisandro Dalcin 453b07a2398SLisandro Dalcin PetscFunctionBegin; 454b07a2398SLisandro Dalcin if (radius < 0 || radius > 1) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Radius %g not in range [0,1]",(double)radius); 455b07a2398SLisandro Dalcin alpha_m = (PetscReal)0.5*(3-radius)/(1+radius); 456b07a2398SLisandro Dalcin alpha_f = 1/(1+radius); 457b07a2398SLisandro Dalcin gamma = (PetscReal)0.5 + alpha_m - alpha_f; 458b07a2398SLisandro Dalcin ierr = TSAlphaSetParams(ts,alpha_m,alpha_f,gamma);CHKERRQ(ierr); 459b07a2398SLisandro Dalcin PetscFunctionReturn(0); 460b07a2398SLisandro Dalcin } 461b07a2398SLisandro Dalcin 462b07a2398SLisandro Dalcin #undef __FUNCT__ 463b07a2398SLisandro Dalcin #define __FUNCT__ "TSAlphaSetParams_Alpha" 464b07a2398SLisandro Dalcin static PetscErrorCode TSAlphaSetParams_Alpha(TS ts,PetscReal alpha_m,PetscReal alpha_f,PetscReal gamma) 465b07a2398SLisandro Dalcin { 466b07a2398SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 467b07a2398SLisandro Dalcin PetscReal tol = 100*PETSC_MACHINE_EPSILON; 468b07a2398SLisandro Dalcin PetscReal res = ((PetscReal)0.5 + alpha_m - alpha_f) - gamma; 469b07a2398SLisandro Dalcin 470b07a2398SLisandro Dalcin PetscFunctionBegin; 471b07a2398SLisandro Dalcin th->Alpha_m = alpha_m; 472b07a2398SLisandro Dalcin th->Alpha_f = alpha_f; 473b07a2398SLisandro Dalcin th->Gamma = gamma; 474b07a2398SLisandro Dalcin th->order = (PetscAbsReal(res) < tol) ? 2 : 1; 475b07a2398SLisandro Dalcin PetscFunctionReturn(0); 476b07a2398SLisandro Dalcin } 477b07a2398SLisandro Dalcin 478b07a2398SLisandro Dalcin #undef __FUNCT__ 479b07a2398SLisandro Dalcin #define __FUNCT__ "TSAlphaGetParams_Alpha" 480b07a2398SLisandro Dalcin static PetscErrorCode TSAlphaGetParams_Alpha(TS ts,PetscReal *alpha_m,PetscReal *alpha_f,PetscReal *gamma) 481b07a2398SLisandro Dalcin { 482b07a2398SLisandro Dalcin TS_Alpha *th = (TS_Alpha*)ts->data; 483b07a2398SLisandro Dalcin 484b07a2398SLisandro Dalcin PetscFunctionBegin; 485b07a2398SLisandro Dalcin if (alpha_m) *alpha_m = th->Alpha_m; 486b07a2398SLisandro Dalcin if (alpha_f) *alpha_f = th->Alpha_f; 487b07a2398SLisandro Dalcin if (gamma) *gamma = th->Gamma; 488b07a2398SLisandro Dalcin PetscFunctionReturn(0); 489b07a2398SLisandro Dalcin } 490b07a2398SLisandro Dalcin 491b07a2398SLisandro Dalcin /*MC 492b07a2398SLisandro Dalcin TSALPHA - ODE/DAE solver using the implicit Generalized-Alpha method 493b07a2398SLisandro Dalcin for first-order systems 494b07a2398SLisandro Dalcin 495b07a2398SLisandro Dalcin Level: beginner 496b07a2398SLisandro Dalcin 497b07a2398SLisandro Dalcin References: 498b07a2398SLisandro Dalcin K.E. Jansen, C.H. Whiting, G.M. Hulber, "A generalized-alpha 499b07a2398SLisandro Dalcin method for integrating the filtered Navier-Stokes equations with a 500b07a2398SLisandro Dalcin stabilized finite element method", Computer Methods in Applied 501b07a2398SLisandro Dalcin Mechanics and Engineering, 190, 305-319, 2000. 502b07a2398SLisandro Dalcin DOI: 10.1016/S0045-7825(00)00203-6. 503b07a2398SLisandro Dalcin 504b07a2398SLisandro Dalcin J. Chung, G.M.Hubert. "A Time Integration Algorithm for Structural 505b07a2398SLisandro Dalcin Dynamics with Improved Numerical Dissipation: The Generalized-alpha 506b07a2398SLisandro Dalcin Method" ASME Journal of Applied Mechanics, 60, 371:375, 1993. 507b07a2398SLisandro Dalcin 508b07a2398SLisandro Dalcin .seealso: TS, TSCreate(), TSSetType(), TSAlphaSetRadius(), TSAlphaSetParams() 509b07a2398SLisandro Dalcin M*/ 510b07a2398SLisandro Dalcin #undef __FUNCT__ 511b07a2398SLisandro Dalcin #define __FUNCT__ "TSCreate_Alpha" 512b07a2398SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSCreate_Alpha(TS ts) 513b07a2398SLisandro Dalcin { 514b07a2398SLisandro Dalcin TS_Alpha *th; 515b07a2398SLisandro Dalcin PetscErrorCode ierr; 516b07a2398SLisandro Dalcin 517b07a2398SLisandro Dalcin PetscFunctionBegin; 518b07a2398SLisandro Dalcin ts->ops->reset = TSReset_Alpha; 519b07a2398SLisandro Dalcin ts->ops->destroy = TSDestroy_Alpha; 520b07a2398SLisandro Dalcin ts->ops->view = TSView_Alpha; 521b07a2398SLisandro Dalcin ts->ops->setup = TSSetUp_Alpha; 522b07a2398SLisandro Dalcin ts->ops->setfromoptions = TSSetFromOptions_Alpha; 523b07a2398SLisandro Dalcin ts->ops->step = TSStep_Alpha; 5249808bdc1SLisandro Dalcin ts->ops->evaluatewlte = TSEvaluateWLTE_Alpha; 525b07a2398SLisandro Dalcin ts->ops->rollback = TSRollBack_Alpha; 526b07a2398SLisandro Dalcin ts->ops->interpolate = TSInterpolate_Alpha; 527b07a2398SLisandro Dalcin ts->ops->snesfunction = SNESTSFormFunction_Alpha; 528b07a2398SLisandro Dalcin ts->ops->snesjacobian = SNESTSFormJacobian_Alpha; 529b07a2398SLisandro Dalcin 530b07a2398SLisandro Dalcin ierr = PetscNewLog(ts,&th);CHKERRQ(ierr); 531b07a2398SLisandro Dalcin ts->data = (void*)th; 532b07a2398SLisandro Dalcin 533b07a2398SLisandro Dalcin th->Alpha_m = 0.5; 534b07a2398SLisandro Dalcin th->Alpha_f = 0.5; 535b07a2398SLisandro Dalcin th->Gamma = 0.5; 536b07a2398SLisandro Dalcin th->order = 2; 537b07a2398SLisandro Dalcin 538b07a2398SLisandro Dalcin th->adapt = PETSC_FALSE; 539b07a2398SLisandro Dalcin 540b07a2398SLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)ts,"TSAlphaUseAdapt_C",TSAlphaUseAdapt_Alpha);CHKERRQ(ierr); 541b07a2398SLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)ts,"TSAlphaSetRadius_C",TSAlphaSetRadius_Alpha);CHKERRQ(ierr); 542b07a2398SLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)ts,"TSAlphaSetParams_C",TSAlphaSetParams_Alpha);CHKERRQ(ierr); 543b07a2398SLisandro Dalcin ierr = PetscObjectComposeFunction((PetscObject)ts,"TSAlphaGetParams_C",TSAlphaGetParams_Alpha);CHKERRQ(ierr); 544b07a2398SLisandro Dalcin PetscFunctionReturn(0); 545b07a2398SLisandro Dalcin } 546b07a2398SLisandro Dalcin 547b07a2398SLisandro Dalcin #undef __FUNCT__ 548b07a2398SLisandro Dalcin #define __FUNCT__ "TSAlphaUseAdapt" 549b07a2398SLisandro Dalcin /*@ 550b07a2398SLisandro Dalcin TSAlphaUseAdapt - Use time-step adaptivity with the Alpha method 551b07a2398SLisandro Dalcin 552b07a2398SLisandro Dalcin Logically Collective on TS 553b07a2398SLisandro Dalcin 554b07a2398SLisandro Dalcin Input Parameter: 555b07a2398SLisandro Dalcin + ts - timestepping context 556b07a2398SLisandro Dalcin - use - flag to use adaptivity 557b07a2398SLisandro Dalcin 558b07a2398SLisandro Dalcin Options Database: 559b07a2398SLisandro Dalcin . -ts_alpha_adapt 560b07a2398SLisandro Dalcin 561b07a2398SLisandro Dalcin Level: intermediate 562b07a2398SLisandro Dalcin 563b07a2398SLisandro Dalcin .seealso: TSAdapt, TSADAPTBASIC 564b07a2398SLisandro Dalcin @*/ 565b07a2398SLisandro Dalcin PetscErrorCode TSAlphaUseAdapt(TS ts,PetscBool use) 566b07a2398SLisandro Dalcin { 567b07a2398SLisandro Dalcin PetscErrorCode ierr; 568b07a2398SLisandro Dalcin 569b07a2398SLisandro Dalcin PetscFunctionBegin; 570b07a2398SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 571b07a2398SLisandro Dalcin PetscValidLogicalCollectiveBool(ts,use,2); 572b07a2398SLisandro Dalcin ierr = PetscTryMethod(ts,"TSAlphaUseAdapt_C",(TS,PetscBool),(ts,use));CHKERRQ(ierr); 573b07a2398SLisandro Dalcin PetscFunctionReturn(0); 574b07a2398SLisandro Dalcin } 575b07a2398SLisandro Dalcin 576b07a2398SLisandro Dalcin #undef __FUNCT__ 577b07a2398SLisandro Dalcin #define __FUNCT__ "TSAlphaSetRadius" 578b07a2398SLisandro Dalcin /*@ 579b07a2398SLisandro Dalcin TSAlphaSetRadius - sets the desired spectral radius of the method 580b07a2398SLisandro Dalcin (i.e. high-frequency numerical damping) 581b07a2398SLisandro Dalcin 582b07a2398SLisandro Dalcin Logically Collective on TS 583b07a2398SLisandro Dalcin 584b07a2398SLisandro Dalcin The algorithmic parameters \alpha_m and \alpha_f of the 585b07a2398SLisandro Dalcin generalized-\alpha method can be computed in terms of a specified 586b07a2398SLisandro Dalcin spectral radius \rho in [0,1] for infinite time step in order to 587b07a2398SLisandro Dalcin control high-frequency numerical damping: 588b07a2398SLisandro Dalcin \alpha_m = 0.5*(3-\rho)/(1+\rho) 589b07a2398SLisandro Dalcin \alpha_f = 1/(1+\rho) 590b07a2398SLisandro Dalcin 591b07a2398SLisandro Dalcin Input Parameter: 592b07a2398SLisandro Dalcin + ts - timestepping context 593b07a2398SLisandro Dalcin - radius - the desired spectral radius 594b07a2398SLisandro Dalcin 595b07a2398SLisandro Dalcin Options Database: 596b07a2398SLisandro Dalcin . -ts_alpha_radius <radius> 597b07a2398SLisandro Dalcin 598b07a2398SLisandro Dalcin Level: intermediate 599b07a2398SLisandro Dalcin 600b07a2398SLisandro Dalcin .seealso: TSAlphaSetParams(), TSAlphaGetParams() 601b07a2398SLisandro Dalcin @*/ 602b07a2398SLisandro Dalcin PetscErrorCode TSAlphaSetRadius(TS ts,PetscReal radius) 603b07a2398SLisandro Dalcin { 604b07a2398SLisandro Dalcin PetscErrorCode ierr; 605b07a2398SLisandro Dalcin 606b07a2398SLisandro Dalcin PetscFunctionBegin; 607b07a2398SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 608b07a2398SLisandro Dalcin PetscValidLogicalCollectiveReal(ts,radius,2); 609b07a2398SLisandro Dalcin if (radius < 0 || radius > 1) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Radius %g not in range [0,1]",(double)radius); 610b07a2398SLisandro Dalcin ierr = PetscTryMethod(ts,"TSAlphaSetRadius_C",(TS,PetscReal),(ts,radius));CHKERRQ(ierr); 611b07a2398SLisandro Dalcin PetscFunctionReturn(0); 612b07a2398SLisandro Dalcin } 613b07a2398SLisandro Dalcin 614b07a2398SLisandro Dalcin #undef __FUNCT__ 615b07a2398SLisandro Dalcin #define __FUNCT__ "TSAlphaSetParams" 616b07a2398SLisandro Dalcin /*@ 617b07a2398SLisandro Dalcin TSAlphaSetParams - sets the algorithmic parameters for TSALPHA 618b07a2398SLisandro Dalcin 619b07a2398SLisandro Dalcin Logically Collective on TS 620b07a2398SLisandro Dalcin 621b07a2398SLisandro Dalcin Second-order accuracy can be obtained so long as: 622b07a2398SLisandro Dalcin \gamma = 0.5 + alpha_m - alpha_f 623b07a2398SLisandro Dalcin 624b07a2398SLisandro Dalcin Unconditional stability requires: 625b07a2398SLisandro Dalcin \alpha_m >= \alpha_f >= 0.5 626b07a2398SLisandro Dalcin 627b07a2398SLisandro Dalcin Backward Euler method is recovered with: 628b07a2398SLisandro Dalcin \alpha_m = \alpha_f = gamma = 1 629b07a2398SLisandro Dalcin 630b07a2398SLisandro Dalcin Input Parameter: 631b07a2398SLisandro Dalcin + ts - timestepping context 632b07a2398SLisandro Dalcin . \alpha_m - algorithmic paramenter 633b07a2398SLisandro Dalcin . \alpha_f - algorithmic paramenter 634b07a2398SLisandro Dalcin - \gamma - algorithmic paramenter 635b07a2398SLisandro Dalcin 636b07a2398SLisandro Dalcin Options Database: 637b07a2398SLisandro Dalcin + -ts_alpha_alpha_m <alpha_m> 638b07a2398SLisandro Dalcin . -ts_alpha_alpha_f <alpha_f> 639b07a2398SLisandro Dalcin - -ts_alpha_gamma <gamma> 640b07a2398SLisandro Dalcin 641b07a2398SLisandro Dalcin Note: 642b07a2398SLisandro Dalcin Use of this function is normally only required to hack TSALPHA to 643b07a2398SLisandro Dalcin use a modified integration scheme. Users should call 644b07a2398SLisandro Dalcin TSAlphaSetRadius() to set the desired spectral radius of the methods 645b07a2398SLisandro Dalcin (i.e. high-frequency damping) in order so select optimal values for 646b07a2398SLisandro Dalcin these parameters. 647b07a2398SLisandro Dalcin 648b07a2398SLisandro Dalcin Level: advanced 649b07a2398SLisandro Dalcin 650b07a2398SLisandro Dalcin .seealso: TSAlphaSetRadius(), TSAlphaGetParams() 651b07a2398SLisandro Dalcin @*/ 652b07a2398SLisandro Dalcin PetscErrorCode TSAlphaSetParams(TS ts,PetscReal alpha_m,PetscReal alpha_f,PetscReal gamma) 653b07a2398SLisandro Dalcin { 654b07a2398SLisandro Dalcin PetscErrorCode ierr; 655b07a2398SLisandro Dalcin 656b07a2398SLisandro Dalcin PetscFunctionBegin; 657b07a2398SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 658b07a2398SLisandro Dalcin PetscValidLogicalCollectiveReal(ts,alpha_m,2); 659b07a2398SLisandro Dalcin PetscValidLogicalCollectiveReal(ts,alpha_f,3); 660b07a2398SLisandro Dalcin PetscValidLogicalCollectiveReal(ts,gamma,4); 661b07a2398SLisandro Dalcin ierr = PetscTryMethod(ts,"TSAlphaSetParams_C",(TS,PetscReal,PetscReal,PetscReal),(ts,alpha_m,alpha_f,gamma));CHKERRQ(ierr); 662b07a2398SLisandro Dalcin PetscFunctionReturn(0); 663b07a2398SLisandro Dalcin } 664b07a2398SLisandro Dalcin 665b07a2398SLisandro Dalcin #undef __FUNCT__ 666b07a2398SLisandro Dalcin #define __FUNCT__ "TSAlphaGetParams" 667b07a2398SLisandro Dalcin /*@ 668b07a2398SLisandro Dalcin TSAlphaGetParams - gets the algorithmic parameters for TSALPHA 669b07a2398SLisandro Dalcin 670b07a2398SLisandro Dalcin Not Collective 671b07a2398SLisandro Dalcin 672b07a2398SLisandro Dalcin Input Parameter: 673b07a2398SLisandro Dalcin . ts - timestepping context 674b07a2398SLisandro Dalcin 675b07a2398SLisandro Dalcin Output Parameters: 676b07a2398SLisandro Dalcin + \alpha_m - algorithmic parameter 677b07a2398SLisandro Dalcin . \alpha_f - algorithmic parameter 678b07a2398SLisandro Dalcin - \gamma - algorithmic parameter 679b07a2398SLisandro Dalcin 680b07a2398SLisandro Dalcin Note: 681b07a2398SLisandro Dalcin Use of this function is normally only required to hack TSALPHA to 682b07a2398SLisandro Dalcin use a modified integration scheme. Users should call 683b07a2398SLisandro Dalcin TSAlphaSetRadius() to set the high-frequency damping (i.e. spectral 684b07a2398SLisandro Dalcin radius of the method) in order so select optimal values for these 685b07a2398SLisandro Dalcin parameters. 686b07a2398SLisandro Dalcin 687b07a2398SLisandro Dalcin Level: advanced 688b07a2398SLisandro Dalcin 689b07a2398SLisandro Dalcin .seealso: TSAlphaSetRadius(), TSAlphaSetParams() 690b07a2398SLisandro Dalcin @*/ 691b07a2398SLisandro Dalcin PetscErrorCode TSAlphaGetParams(TS ts,PetscReal *alpha_m,PetscReal *alpha_f,PetscReal *gamma) 692b07a2398SLisandro Dalcin { 693b07a2398SLisandro Dalcin PetscErrorCode ierr; 694b07a2398SLisandro Dalcin 695b07a2398SLisandro Dalcin PetscFunctionBegin; 696b07a2398SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,1); 697b07a2398SLisandro Dalcin if (alpha_m) PetscValidRealPointer(alpha_m,2); 698b07a2398SLisandro Dalcin if (alpha_f) PetscValidRealPointer(alpha_f,3); 699b07a2398SLisandro Dalcin if (gamma) PetscValidRealPointer(gamma,4); 700b07a2398SLisandro Dalcin ierr = PetscUseMethod(ts,"TSAlphaGetParams_C",(TS,PetscReal*,PetscReal*,PetscReal*),(ts,alpha_m,alpha_f,gamma));CHKERRQ(ierr); 701b07a2398SLisandro Dalcin PetscFunctionReturn(0); 702b07a2398SLisandro Dalcin } 703