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