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