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