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