Actual source code: alpha1.c
petsc-3.7.3 2016-08-01
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*/
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";
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;
27: PetscReal Alpha_m;
28: PetscReal Alpha_f;
29: PetscReal Gamma;
30: PetscInt order;
32: PetscBool adapt;
33: Vec vec_sol_prev;
34: Vec vec_lte_work;
36: TSStepStatus status;
37: } TS_Alpha;
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;
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: return(0);
55: }
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;
72: /* V1 = 1/(Gamma*dT)*(X1-X0) + (1-1/Gamma)*V0 */
73: VecWAXPY(V1,-1.0,X0,X1);
74: VecAXPBY(V1,1-1/Gamma,1/(Gamma*dt),V0);
75: /* Xa = X0 + Alpha_f*(X1-X0) */
76: VecWAXPY(Xa,-1.0,X0,X1);
77: VecAYPX(Xa,Alpha_f,X0);
78: /* Va = V0 + Alpha_m*(V1-V0) */
79: VecWAXPY(Va,-1.0,V0,V1);
80: VecAYPX(Va,Alpha_m,V0);
81: return(0);
82: }
86: static PetscErrorCode TS_SNESSolve(TS ts,Vec b,Vec x)
87: {
88: PetscInt nits,lits;
92: SNESSolve(ts->snes,b,x);
93: SNESGetIterationNumber(ts->snes,&nits);
94: SNESGetLinearSolveIterations(ts->snes,&lits);
95: ts->snes_its += nits; ts->ksp_its += lits;
96: return(0);
97: }
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: */
107: static PetscErrorCode TSAlpha_Restart(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;
117: VecDuplicate(X0,&X1);
119: /* Setup backward Euler with halved time step */
120: TSAlphaGetParams(ts,&alpha_m,&alpha_f,&gamma);
121: TSAlphaSetParams(ts,1,1,1);
122: TSGetTimeStep(ts,&time_step);
123: ts->time_step = time_step/2;
124: TSAlpha_StageTime(ts);
125: th->stage_time = ts->ptime;
126: VecZeroEntries(th->V0);
128: /* First BE step, (t0,X0) -> (t1,X1) */
129: th->stage_time += ts->time_step;
130: VecCopy(X0,th->X0);
131: TSPreStage(ts,th->stage_time);
132: VecCopy(th->X0,X1);
133: TS_SNESSolve(ts,NULL,X1);
134: TSPostStage(ts,th->stage_time,0,&X1);
135: TSAdaptCheckStage(ts->adapt,ts,th->stage_time,X1,&stageok);
136: if (!stageok) goto finally;
138: /* Second BE step, (t1,X1) -> (t2,X2) */
139: th->stage_time += ts->time_step;
140: VecCopy(X1,th->X0);
141: TSPreStage(ts,th->stage_time);
142: VecCopy(th->X0,X2);
143: TS_SNESSolve(ts,NULL,X2);
144: TSPostStage(ts,th->stage_time,0,&X2);
145: TSAdaptCheckStage(ts->adapt,ts,th->stage_time,X2,&stageok);
146: if (!stageok) goto finally;
148: /* Compute V0 ~ dX/dt at t0 with backward differences */
149: VecZeroEntries(th->V0);
150: VecAXPY(th->V0,-3/ts->time_step,X0);
151: VecAXPY(th->V0,+4/ts->time_step,X1);
152: VecAXPY(th->V0,-1/ts->time_step,X2);
154: /* Rough, lower-order estimate LTE of the initial step */
155: if (th->adapt) {
156: VecZeroEntries(th->vec_lte_work);
157: VecAXPY(th->vec_lte_work,+2,X2);
158: VecAXPY(th->vec_lte_work,-4,X1);
159: VecAXPY(th->vec_lte_work,+2,X0);
160: }
162: finally:
163: /* Revert TSAlpha to the initial state (t0,X0) */
164: if (initok) *initok = stageok;
165: TSSetTimeStep(ts,time_step);
166: TSAlphaSetParams(ts,alpha_m,alpha_f,gamma);
167: VecCopy(ts->vec_sol,th->X0);
169: VecDestroy(&X1);
170: return(0);
171: }
175: static PetscErrorCode TSStep_Alpha(TS ts)
176: {
177: TS_Alpha *th = (TS_Alpha*)ts->data;
178: PetscInt rejections = 0;
179: PetscBool stageok,accept = PETSC_TRUE;
180: PetscReal next_time_step = ts->time_step;
184: PetscCitationsRegister(citation,&cited);
186: if (!ts->steprollback) {
187: if (th->adapt) { VecCopy(th->X0,th->vec_sol_prev); }
188: VecCopy(ts->vec_sol,th->X0);
189: VecCopy(th->V1,th->V0);
190: }
192: th->status = TS_STEP_INCOMPLETE;
193: while (!ts->reason && th->status != TS_STEP_COMPLETE) {
195: if (ts->steprestart) {
196: TSAlpha_Restart(ts,&stageok);
197: if (!stageok) goto reject_step;
198: }
200: TSAlpha_StageTime(ts);
201: VecCopy(th->X0,th->X1);
202: TSPreStage(ts,th->stage_time);
203: TS_SNESSolve(ts,NULL,th->X1);
204: TSPostStage(ts,th->stage_time,0,&th->Xa);
205: TSAdaptCheckStage(ts->adapt,ts,th->stage_time,th->Xa,&stageok);
206: if (!stageok) goto reject_step;
208: th->status = TS_STEP_PENDING;
209: VecCopy(th->X1,ts->vec_sol);
210: TSAdaptChoose(ts->adapt,ts,ts->time_step,NULL,&next_time_step,&accept);
211: th->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
212: if (!accept) {
213: VecCopy(th->X0,ts->vec_sol);
214: ts->time_step = next_time_step;
215: goto reject_step;
216: }
218: ts->ptime += ts->time_step;
219: ts->time_step = next_time_step;
220: break;
222: reject_step:
223: ts->reject++; accept = PETSC_FALSE;
224: if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
225: ts->reason = TS_DIVERGED_STEP_REJECTED;
226: PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);
227: }
229: }
230: return(0);
231: }
235: static PetscErrorCode TSEvaluateWLTE_Alpha(TS ts,NormType wnormtype,PetscInt *order,PetscReal *wlte)
236: {
237: TS_Alpha *th = (TS_Alpha*)ts->data;
238: Vec X = th->X1; /* X = solution */
239: Vec Y = th->vec_lte_work; /* Y = X + LTE */
243: if (ts->steprestart) {
244: /* th->vec_lte_work is set to the LTE in TSAlpha_Restart() */
245: VecAXPY(Y,1,X);
246: } else {
247: /* Compute LTE using backward differences with non-constant time step */
248: PetscReal h = ts->time_step, h_prev = ts->ptime - ts->ptime_prev;
249: PetscReal a = 1 + h_prev/h;
250: PetscScalar scal[3]; Vec vecs[3];
251: scal[0] = +1/a; scal[1] = -1/(a-1); scal[2] = +1/(a*(a-1));
252: vecs[0] = th->X1; vecs[1] = th->X0; vecs[2] = th->vec_sol_prev;
253: VecCopy(X,Y);
254: VecMAXPY(Y,3,scal,vecs);
255: }
256: TSErrorWeightedNorm(ts,X,Y,wnormtype,wlte);
257: if (order) *order = 2;
258: return(0);
259: }
263: static PetscErrorCode TSRollBack_Alpha(TS ts)
264: {
265: TS_Alpha *th = (TS_Alpha*)ts->data;
269: VecCopy(th->X0,ts->vec_sol);
270: return(0);
271: }
275: static PetscErrorCode TSInterpolate_Alpha(TS ts,PetscReal t,Vec X)
276: {
277: TS_Alpha *th = (TS_Alpha*)ts->data;
278: PetscReal dt = t - ts->ptime;
282: VecCopy(ts->vec_sol,X);
283: VecAXPY(X,th->Gamma*dt,th->V1);
284: VecAXPY(X,(1-th->Gamma)*dt,th->V0);
285: return(0);
286: }
290: static PetscErrorCode SNESTSFormFunction_Alpha(PETSC_UNUSED SNES snes,Vec X,Vec F,TS ts)
291: {
292: TS_Alpha *th = (TS_Alpha*)ts->data;
293: PetscReal ta = th->stage_time;
294: Vec Xa = th->Xa, Va = th->Va;
298: TSAlpha_StageVecs(ts,X);
299: /* F = Function(ta,Xa,Va) */
300: TSComputeIFunction(ts,ta,Xa,Va,F,PETSC_FALSE);
301: VecScale(F,th->scale_F);
302: return(0);
303: }
307: static PetscErrorCode SNESTSFormJacobian_Alpha(PETSC_UNUSED SNES snes,PETSC_UNUSED Vec X,Mat J,Mat P,TS ts)
308: {
309: TS_Alpha *th = (TS_Alpha*)ts->data;
310: PetscReal ta = th->stage_time;
311: Vec Xa = th->Xa, Va = th->Va;
312: PetscReal dVdX = th->shift_V;
316: /* J,P = Jacobian(ta,Xa,Va) */
317: TSComputeIJacobian(ts,ta,Xa,Va,dVdX,J,P,PETSC_FALSE);
318: return(0);
319: }
323: static PetscErrorCode TSReset_Alpha(TS ts)
324: {
325: TS_Alpha *th = (TS_Alpha*)ts->data;
329: VecDestroy(&th->X0);
330: VecDestroy(&th->Xa);
331: VecDestroy(&th->X1);
332: VecDestroy(&th->V0);
333: VecDestroy(&th->Va);
334: VecDestroy(&th->V1);
335: VecDestroy(&th->vec_sol_prev);
336: VecDestroy(&th->vec_lte_work);
337: return(0);
338: }
342: static PetscErrorCode TSDestroy_Alpha(TS ts)
343: {
347: TSReset_Alpha(ts);
348: PetscFree(ts->data);
350: PetscObjectComposeFunction((PetscObject)ts,"TSAlphaUseAdapt_C",NULL);
351: PetscObjectComposeFunction((PetscObject)ts,"TSAlphaSetRadius_C",NULL);
352: PetscObjectComposeFunction((PetscObject)ts,"TSAlphaSetParams_C",NULL);
353: PetscObjectComposeFunction((PetscObject)ts,"TSAlphaGetParams_C",NULL);
354: return(0);
355: }
359: static PetscErrorCode TSSetUp_Alpha(TS ts)
360: {
361: TS_Alpha *th = (TS_Alpha*)ts->data;
365: VecDuplicate(ts->vec_sol,&th->X0);
366: VecDuplicate(ts->vec_sol,&th->Xa);
367: VecDuplicate(ts->vec_sol,&th->X1);
368: VecDuplicate(ts->vec_sol,&th->V0);
369: VecDuplicate(ts->vec_sol,&th->Va);
370: VecDuplicate(ts->vec_sol,&th->V1);
372: TSGetAdapt(ts,&ts->adapt);
373: TSAdaptCandidatesClear(ts->adapt);
374: if (!th->adapt) {
375: TSAdaptSetType(ts->adapt,TSADAPTNONE);
376: } else {
377: VecDuplicate(ts->vec_sol,&th->vec_sol_prev);
378: VecDuplicate(ts->vec_sol,&th->vec_lte_work);
379: if (ts->exact_final_time == TS_EXACTFINALTIME_UNSPECIFIED)
380: ts->exact_final_time = TS_EXACTFINALTIME_MATCHSTEP;
381: }
383: TSGetSNES(ts,&ts->snes);
384: return(0);
385: }
389: static PetscErrorCode TSSetFromOptions_Alpha(PetscOptionItems *PetscOptionsObject,TS ts)
390: {
391: TS_Alpha *th = (TS_Alpha*)ts->data;
395: PetscOptionsHead(PetscOptionsObject,"Generalized-Alpha ODE solver options");
396: {
397: PetscBool flg;
398: PetscReal radius = 1;
399: PetscBool adapt = th->adapt;
400: PetscOptionsReal("-ts_alpha_radius","Spectral radius (high-frequency dissipation)","TSAlphaSetRadius",radius,&radius,&flg);
401: if (flg) {TSAlphaSetRadius(ts,radius);}
402: PetscOptionsReal("-ts_alpha_alpha_m","Algoritmic parameter alpha_m","TSAlphaSetParams",th->Alpha_m,&th->Alpha_m,NULL);
403: PetscOptionsReal("-ts_alpha_alpha_f","Algoritmic parameter alpha_f","TSAlphaSetParams",th->Alpha_f,&th->Alpha_f,NULL);
404: PetscOptionsReal("-ts_alpha_gamma","Algoritmic parameter gamma","TSAlphaSetParams",th->Gamma,&th->Gamma,NULL);
405: TSAlphaSetParams(ts,th->Alpha_m,th->Alpha_f,th->Gamma);
406: PetscOptionsBool("-ts_alpha_adapt","Use time-step adaptivity with the Alpha method","TSAlpha2UseAdapt",adapt,&adapt,&flg);
407: if (flg) {TSAlphaUseAdapt(ts,adapt);}
408: }
409: PetscOptionsTail();
410: return(0);
411: }
415: static PetscErrorCode TSView_Alpha(TS ts,PetscViewer viewer)
416: {
417: TS_Alpha *th = (TS_Alpha*)ts->data;
418: PetscBool iascii;
422: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
423: if (iascii) {
424: PetscViewerASCIIPrintf(viewer," Alpha_m=%g, Alpha_f=%g, Gamma=%g\n",(double)th->Alpha_m,(double)th->Alpha_f,(double)th->Gamma);
425: }
426: if (ts->adapt) {TSAdaptView(ts->adapt,viewer);}
427: if (ts->snes) {SNESView(ts->snes,viewer);}
428: return(0);
429: }
433: static PetscErrorCode TSAlphaUseAdapt_Alpha(TS ts,PetscBool use)
434: {
435: TS_Alpha *th = (TS_Alpha*)ts->data;
438: if (use == th->adapt) return(0);
439: if (ts->setupcalled) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ORDER,"Cannot change adaptivity after TSSetUp()");
440: th->adapt = use;
441: return(0);
442: }
446: static PetscErrorCode TSAlphaSetRadius_Alpha(TS ts,PetscReal radius)
447: {
448: PetscReal alpha_m,alpha_f,gamma;
452: if (radius < 0 || radius > 1) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Radius %g not in range [0,1]",(double)radius);
453: alpha_m = (PetscReal)0.5*(3-radius)/(1+radius);
454: alpha_f = 1/(1+radius);
455: gamma = (PetscReal)0.5 + alpha_m - alpha_f;
456: TSAlphaSetParams(ts,alpha_m,alpha_f,gamma);
457: return(0);
458: }
462: static PetscErrorCode TSAlphaSetParams_Alpha(TS ts,PetscReal alpha_m,PetscReal alpha_f,PetscReal gamma)
463: {
464: TS_Alpha *th = (TS_Alpha*)ts->data;
465: PetscReal tol = 100*PETSC_MACHINE_EPSILON;
466: PetscReal res = ((PetscReal)0.5 + alpha_m - alpha_f) - gamma;
469: th->Alpha_m = alpha_m;
470: th->Alpha_f = alpha_f;
471: th->Gamma = gamma;
472: th->order = (PetscAbsReal(res) < tol) ? 2 : 1;
473: return(0);
474: }
478: static PetscErrorCode TSAlphaGetParams_Alpha(TS ts,PetscReal *alpha_m,PetscReal *alpha_f,PetscReal *gamma)
479: {
480: TS_Alpha *th = (TS_Alpha*)ts->data;
483: if (alpha_m) *alpha_m = th->Alpha_m;
484: if (alpha_f) *alpha_f = th->Alpha_f;
485: if (gamma) *gamma = th->Gamma;
486: return(0);
487: }
489: /*MC
490: TSALPHA - ODE/DAE solver using the implicit Generalized-Alpha method
491: for first-order systems
493: Level: beginner
495: References:
496: K.E. Jansen, C.H. Whiting, G.M. Hulber, "A generalized-alpha
497: method for integrating the filtered Navier-Stokes equations with a
498: stabilized finite element method", Computer Methods in Applied
499: Mechanics and Engineering, 190, 305-319, 2000.
500: DOI: 10.1016/S0045-7825(00)00203-6.
502: J. Chung, G.M.Hubert. "A Time Integration Algorithm for Structural
503: Dynamics with Improved Numerical Dissipation: The Generalized-alpha
504: Method" ASME Journal of Applied Mechanics, 60, 371:375, 1993.
506: .seealso: TS, TSCreate(), TSSetType(), TSAlphaSetRadius(), TSAlphaSetParams()
507: M*/
510: PETSC_EXTERN PetscErrorCode TSCreate_Alpha(TS ts)
511: {
512: TS_Alpha *th;
516: ts->ops->reset = TSReset_Alpha;
517: ts->ops->destroy = TSDestroy_Alpha;
518: ts->ops->view = TSView_Alpha;
519: ts->ops->setup = TSSetUp_Alpha;
520: ts->ops->setfromoptions = TSSetFromOptions_Alpha;
521: ts->ops->step = TSStep_Alpha;
522: ts->ops->evaluatewlte = TSEvaluateWLTE_Alpha;
523: ts->ops->rollback = TSRollBack_Alpha;
524: ts->ops->interpolate = TSInterpolate_Alpha;
525: ts->ops->snesfunction = SNESTSFormFunction_Alpha;
526: ts->ops->snesjacobian = SNESTSFormJacobian_Alpha;
528: PetscNewLog(ts,&th);
529: ts->data = (void*)th;
531: th->Alpha_m = 0.5;
532: th->Alpha_f = 0.5;
533: th->Gamma = 0.5;
534: th->order = 2;
536: th->adapt = PETSC_FALSE;
538: PetscObjectComposeFunction((PetscObject)ts,"TSAlphaUseAdapt_C",TSAlphaUseAdapt_Alpha);
539: PetscObjectComposeFunction((PetscObject)ts,"TSAlphaSetRadius_C",TSAlphaSetRadius_Alpha);
540: PetscObjectComposeFunction((PetscObject)ts,"TSAlphaSetParams_C",TSAlphaSetParams_Alpha);
541: PetscObjectComposeFunction((PetscObject)ts,"TSAlphaGetParams_C",TSAlphaGetParams_Alpha);
542: return(0);
543: }
547: /*@
548: TSAlphaUseAdapt - Use time-step adaptivity with the Alpha method
550: Logically Collective on TS
552: Input Parameter:
553: + ts - timestepping context
554: - use - flag to use adaptivity
556: Options Database:
557: . -ts_alpha_adapt
559: Level: intermediate
561: .seealso: TSAdapt, TSADAPTBASIC
562: @*/
563: PetscErrorCode TSAlphaUseAdapt(TS ts,PetscBool use)
564: {
570: PetscTryMethod(ts,"TSAlphaUseAdapt_C",(TS,PetscBool),(ts,use));
571: return(0);
572: }
576: /*@
577: TSAlphaSetRadius - sets the desired spectral radius of the method
578: (i.e. high-frequency numerical damping)
580: Logically Collective on TS
582: The algorithmic parameters \alpha_m and \alpha_f of the
583: generalized-\alpha method can be computed in terms of a specified
584: spectral radius \rho in [0,1] for infinite time step in order to
585: control high-frequency numerical damping:
586: \alpha_m = 0.5*(3-\rho)/(1+\rho)
587: \alpha_f = 1/(1+\rho)
589: Input Parameter:
590: + ts - timestepping context
591: - radius - the desired spectral radius
593: Options Database:
594: . -ts_alpha_radius <radius>
596: Level: intermediate
598: .seealso: TSAlphaSetParams(), TSAlphaGetParams()
599: @*/
600: PetscErrorCode TSAlphaSetRadius(TS ts,PetscReal radius)
601: {
607: if (radius < 0 || radius > 1) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Radius %g not in range [0,1]",(double)radius);
608: PetscTryMethod(ts,"TSAlphaSetRadius_C",(TS,PetscReal),(ts,radius));
609: return(0);
610: }
614: /*@
615: TSAlphaSetParams - sets the algorithmic parameters for TSALPHA
617: Logically Collective on TS
619: Second-order accuracy can be obtained so long as:
620: \gamma = 0.5 + alpha_m - alpha_f
622: Unconditional stability requires:
623: \alpha_m >= \alpha_f >= 0.5
625: Backward Euler method is recovered with:
626: \alpha_m = \alpha_f = gamma = 1
628: Input Parameter:
629: + ts - timestepping context
630: . \alpha_m - algorithmic paramenter
631: . \alpha_f - algorithmic paramenter
632: - \gamma - algorithmic paramenter
634: Options Database:
635: + -ts_alpha_alpha_m <alpha_m>
636: . -ts_alpha_alpha_f <alpha_f>
637: - -ts_alpha_gamma <gamma>
639: Note:
640: Use of this function is normally only required to hack TSALPHA to
641: use a modified integration scheme. Users should call
642: TSAlphaSetRadius() to set the desired spectral radius of the methods
643: (i.e. high-frequency damping) in order so select optimal values for
644: these parameters.
646: Level: advanced
648: .seealso: TSAlphaSetRadius(), TSAlphaGetParams()
649: @*/
650: PetscErrorCode TSAlphaSetParams(TS ts,PetscReal alpha_m,PetscReal alpha_f,PetscReal gamma)
651: {
659: PetscTryMethod(ts,"TSAlphaSetParams_C",(TS,PetscReal,PetscReal,PetscReal),(ts,alpha_m,alpha_f,gamma));
660: return(0);
661: }
665: /*@
666: TSAlphaGetParams - gets the algorithmic parameters for TSALPHA
668: Not Collective
670: Input Parameter:
671: . ts - timestepping context
673: Output Parameters:
674: + \alpha_m - algorithmic parameter
675: . \alpha_f - algorithmic parameter
676: - \gamma - algorithmic parameter
678: Note:
679: Use of this function is normally only required to hack TSALPHA to
680: use a modified integration scheme. Users should call
681: TSAlphaSetRadius() to set the high-frequency damping (i.e. spectral
682: radius of the method) in order so select optimal values for these
683: parameters.
685: Level: advanced
687: .seealso: TSAlphaSetRadius(), TSAlphaSetParams()
688: @*/
689: PetscErrorCode TSAlphaGetParams(TS ts,PetscReal *alpha_m,PetscReal *alpha_f,PetscReal *gamma)
690: {
698: PetscUseMethod(ts,"TSAlphaGetParams_C",(TS,PetscReal*,PetscReal*,PetscReal*),(ts,alpha_m,alpha_f,gamma));
699: return(0);
700: }