4/5th Order Dormand–Prince Runge-Kutta ODE solver. More...
Public Member Functions | |
TypeName ("RKDP45") | |
Runtime type information. More... | |
RKDP45 (const ODESystem &ode, const dictionary &dict) | |
Construct from ODE. More... | |
scalar | solve (const scalar x0, const scalarField &y0, const scalarField &dydx0, const scalar dx, scalarField &y) const |
Solve a single step dx and return the error. More... | |
void | solve (scalar &x, scalarField &y, scalar &dxTry) const |
Solve the ODE system and the update the state. More... | |
virtual void | solve (scalar &x, scalarField &y, scalar &dxTry) const |
Inherit solve from ODESolver. More... | |
virtual void | solve (scalar &x, scalarField &y, stepState &step) const |
Inherit solve from ODESolver. More... | |
virtual void | solve (const scalar xStart, const scalar xEnd, scalarField &y, scalar &dxEst) const |
Inherit solve from ODESolver. More... | |
![]() | |
TypeName ("ODESolver") | |
Runtime type information. More... | |
declareRunTimeSelectionTable (autoPtr, ODESolver, dictionary,(const ODESystem &ode, const dictionary &dict),(ode, dict)) | |
ODESolver (const ODESystem &ode, const dictionary &dict) | |
Construct for given ODESystem. More... | |
ODESolver (const ODESystem &ode, const scalarField &absTol, const scalarField &relTol) | |
Construct for given ODESystem specifying tolerances. More... | |
virtual | ~ODESolver () |
Destructor. More... | |
scalarField & | absTol () |
scalarField & | relTol () |
virtual void | solve (scalar &x, scalarField &y, stepState &step) const |
Solve the ODE system as far as possible upto dxTry. More... | |
virtual void | solve (const scalar xStart, const scalar xEnd, scalarField &y, scalar &dxEst) const |
Solve the ODE system from xStart to xEnd, update the state. More... | |
![]() | |
adaptiveSolver (const ODESystem &ode, const dictionary &dict) | |
Construct from ODESystem. More... | |
virtual | ~adaptiveSolver () |
Destructor. More... | |
void | solve (const ODESystem &ode, scalar &x, scalarField &y, scalar &dxTry) const |
Solve the ODE system and the update the state. More... | |
Private Attributes | |
scalarField | yTemp_ |
scalarField | k2_ |
scalarField | k3_ |
scalarField | k4_ |
scalarField | k5_ |
scalarField | k6_ |
scalarField | err_ |
Error-estimate field. More... | |
Static Private Attributes | |
static const scalar | c2 = 1.0/5.0 |
RKDP Constants. More... | |
static const scalar | c3 = 3.0/10.0 |
static const scalar | c4 = 4.0/5.0 |
static const scalar | c5 = 8.0/9.0 |
static const scalar | a21 = 1.0/5.0 |
static const scalar | a31 = 3.0/40.0 |
static const scalar | a32 = 9.0/40.0 |
static const scalar | a41 = 44.0/45.0 |
static const scalar | a42 = -56.0/15.0 |
static const scalar | a43 = 32.0/9.0 |
static const scalar | a51 = 19372.0/6561.0 |
static const scalar | a52 = -25360.0/2187.0 |
static const scalar | a53 = 64448.0/6561.0 |
static const scalar | a54 = -212.0/729.0 |
static const scalar | a61 = 9017.0/3168.0 |
static const scalar | a62 = -355.0/33.0 |
static const scalar | a63 = 46732.0/5247.0 |
static const scalar | a64 = 49.0/176.0 |
static const scalar | a65 = -5103.0/18656.0 |
static const scalar | b1 = 35.0/384.0 |
static const scalar | b3 = 500.0/1113.0 |
static const scalar | b4 = 125.0/192.0 |
static const scalar | b5 = -2187.0/6784.0 |
static const scalar | b6 = 11.0/84.0 |
static const scalar | e1 = 5179.0/57600.0 - RKDP45::b1 |
static const scalar | e3 = 7571.0/16695.0 - RKDP45::b3 |
static const scalar | e4 = 393.0/640.0 - RKDP45::b4 |
static const scalar | e5 = -92097.0/339200.0 - RKDP45::b5 |
static const scalar | e6 = 187.0/2100.0 - RKDP45::b6 |
static const scalar | e7 = 1.0/40.0 |
Additional Inherited Members | |
![]() | |
static autoPtr< ODESolver > | New (const ODESystem &ode, const dictionary &dict) |
Select null constructed. More... | |
![]() | |
scalar | normalizeError (const scalarField &y0, const scalarField &y, const scalarField &err) const |
Return the nomalized scalar error. More... | |
ODESolver (const ODESolver &) | |
Disallow default bitwise copy construct. More... | |
void | operator= (const ODESolver &) |
Disallow default bitwise assignment. More... | |
![]() | |
const ODESystem & | odes_ |
Reference to ODESystem. More... | |
label | n_ |
Size of the ODESystem. More... | |
scalarField | absTol_ |
Absolute convergence tolerance per step. More... | |
scalarField | relTol_ |
Relative convergence tolerance per step. More... | |
label | maxSteps_ |
The maximum number of sub-steps allowed for the integration step. More... | |
4/5th Order Dormand–Prince Runge-Kutta ODE solver.
"A family of embedded Runge-Kutta formulae" Dormand, J. R., Prince, P. J., Journal of Computational and Applied Mathematics, 6 (1), 1980: pp. 19-26. "Solving Ordinary Differential Equations I: Nonstiff Problems, second edition", Hairer, E., Nørsett, S., Wanner, G., Springer-Verlag, Berlin. 1993, ISBN 3-540-56670-8.
RKDP45 | ( | const ODESystem & | ode, |
const dictionary & | dict | ||
) |
TypeName | ( | "RKDP45" | ) |
Runtime type information.
|
virtual |
Solve a single step dx and return the error.
Implements adaptiveSolver.
Definition at line 96 of file RKDP45.C.
References Foam::constant::physicoChemical::c2, forAll, y, and Foam::y0().
|
virtual |
|
inline |
Inherit solve from ODESolver.
Definition at line 178 of file ODESolver.H.
void solve |
Inherit solve from ODESolver.
Definition at line 86 of file ODESolver.C.
void solve |
Inherit solve from ODESolver.
Definition at line 99 of file ODESolver.C.
|
staticprivate |
|
staticprivate |
|
staticprivate |
|
staticprivate |
|
staticprivate |
|
mutableprivate |
|
mutableprivate |
|
mutableprivate |
|
mutableprivate |
|
mutableprivate |
|
mutableprivate |
|
mutableprivate |
Copyright © 2011-2018 OpenFOAM | OPENFOAM® is a registered trademark of OpenCFD Ltd.