![]() |
OpenSWMM Engine
6.0.0-alpha.1
Data-oriented, plugin-extensible SWMM Engine (6.0.0-alpha.1)
|
Typedefs | |
| using | DerivFunc = std::function< void(double x, const double *y, double *dydx)> |
| Derivative function signature. | |
| using | BatchDerivFunc = std::function< void(double x, const double *y0, const double *y1, double *dy0, double *dy1, int n_sys)> |
| Integrate N independent 2-equation ODE systems in batch. | |
Functions | |
| int | integrate (double *y, int n, double x1, double x2, double eps, double h1, const DerivFunc &derivs) |
| Integrate an ODE system from x1 to x2 using RK45 Cash-Karp. | |
| int | integrate_batch_2eq (double *y0, double *y1, int n_sys, double x1, double x2, double eps, double h1, const BatchDerivFunc &derivs) |
Variables | |
| constexpr int | MAXSTP = 10000 |
| Maximum integration steps. | |
| constexpr double | TINY = 1.0e-30 |
| Underflow protection. | |
| constexpr double | SAFETY = 0.9 |
| Step adjustment safety factor. | |
| constexpr double | PGROW = -0.2 |
| Exponent for step increase. | |
| constexpr double | PSHRNK = -0.25 |
| Exponent for step decrease. | |
| constexpr double | ERRCON = 1.89e-4 |
| = (5/SAFETY)^(1/PGROW) | |
| constexpr int | ODE_OK = 0 |
| Success. | |
| constexpr int | ODE_TOO_MANY = 1 |
| n > allocated max | |
| constexpr int | ODE_UNDERFLOW = 2 |
| Step size underflowed. | |
| constexpr int | ODE_MAX_STEPS = 3 |
| Exceeded MAXSTP. | |
| using openswmm::ode::BatchDerivFunc = typedef std::function<void(double x, const double* y0, const double* y1, double* dy0, double* dy1, int n_sys)> |
Integrate N independent 2-equation ODE systems in batch.
Used by groundwater: each subcatchment has 2 state variables (theta, lower_depth). All are integrated over the same interval [x1, x2] but with independent derivative functions.
The batch version avoids per-subcatchment function call overhead and enables future SIMD parallelism of the RK stages.
| y0 | [in/out] SoA: first state variable for all N systems. |
| y1 | [in/out] SoA: second state variable for all N systems. |
| n_sys | Number of independent ODE systems. |
| x1 | Start of interval. |
| x2 | End of interval. |
| eps | Tolerance. |
| h1 | Initial step size. |
| derivs | Batch derivative: derivs(x, y0[], y1[], dy0[], dy1[], n_sys). |
| using openswmm::ode::DerivFunc = typedef std::function<void(double x, const double* y, double* dydx)> |
Derivative function signature.
| x | Independent variable (e.g., time). |
| y | Array of n dependent variables at x. |
| dydx | [out] Array where callback writes dy_i/dx. |
| int openswmm::ode::integrate | ( | double * | y, |
| int | n, | ||
| double | x1, | ||
| double | x2, | ||
| double | eps, | ||
| double | h1, | ||
| const DerivFunc & | derivs | ||
| ) |
Integrate an ODE system from x1 to x2 using RK45 Cash-Karp.
| y | [in/out] Array of n dependent variables (initial → final). |
| n | Number of equations. |
| x1 | Start of integration interval. |
| x2 | End of integration interval. |
| eps | Desired accuracy (tolerance). |
| h1 | Initial step size guess. |
| derivs | Derivative callback function. |
| int openswmm::ode::integrate_batch_2eq | ( | double * | y0, |
| double * | y1, | ||
| int | n_sys, | ||
| double | x1, | ||
| double | x2, | ||
| double | eps, | ||
| double | h1, | ||
| const BatchDerivFunc & | derivs | ||
| ) |
|
constexpr |
= (5/SAFETY)^(1/PGROW)
|
constexpr |
Maximum integration steps.
|
constexpr |
Exceeded MAXSTP.
|
constexpr |
Success.
|
constexpr |
n > allocated max
|
constexpr |
Step size underflowed.
|
constexpr |
Exponent for step increase.
|
constexpr |
Exponent for step decrease.
|
constexpr |
Step adjustment safety factor.
|
constexpr |
Underflow protection.