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