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