OpenSWMM Engine  6.0.0-alpha.1
Data-oriented, plugin-extensible SWMM Engine (6.0.0-alpha.1)
Loading...
Searching...
No Matches
OdeSolver.hpp
Go to the documentation of this file.
1
26#ifndef OPENSWMM_ODE_SOLVER_HPP
27#define OPENSWMM_ODE_SOLVER_HPP
28
29#include <functional>
30#include <vector>
31
32namespace openswmm {
33namespace ode {
34
35// ============================================================================
36// Constants (matching legacy)
37// ============================================================================
38
39constexpr int MAXSTP = 10000;
40constexpr double TINY = 1.0e-30;
41constexpr double SAFETY = 0.9;
42constexpr double PGROW = -0.2;
43constexpr double PSHRNK = -0.25;
44constexpr double ERRCON = 1.89e-4;
45
46// ============================================================================
47// Return codes
48// ============================================================================
49
50constexpr int ODE_OK = 0;
51constexpr int ODE_TOO_MANY = 1;
52constexpr int ODE_UNDERFLOW = 2;
53constexpr int ODE_MAX_STEPS = 3;
54
55// ============================================================================
56// Derivative callback
57// ============================================================================
58
66using DerivFunc = std::function<void(double x, const double* y, double* dydx)>;
67
68// ============================================================================
69// Per-element integrator
70// ============================================================================
71
84int integrate(double* y, int n, double x1, double x2,
85 double eps, double h1, const DerivFunc& derivs);
86
87// ============================================================================
88// Batch integrator (SoA, for vectorized subcatchment processing)
89// ============================================================================
90
110using BatchDerivFunc = std::function<void(double x,
111 const double* y0, const double* y1,
112 double* dy0, double* dy1, int n_sys)>;
113
114int integrate_batch_2eq(double* y0, double* y1, int n_sys,
115 double x1, double x2, double eps, double h1,
116 const BatchDerivFunc& derivs);
117
118} // namespace ode
119} // namespace openswmm
120
121#endif // OPENSWMM_ODE_SOLVER_HPP
constexpr int MAXSTP
Maximum integration steps.
Definition OdeSolver.hpp:39
constexpr double PSHRNK
Exponent for step decrease.
Definition OdeSolver.hpp:43
constexpr double SAFETY
Step adjustment safety factor.
Definition OdeSolver.hpp:41
std::function< void(double x, const double *y, double *dydx)> DerivFunc
Derivative function signature.
Definition OdeSolver.hpp:66
constexpr double PGROW
Exponent for step increase.
Definition OdeSolver.hpp:42
std::function< void(double x, const double *y0, const double *y1, double *dy0, double *dy1, int n_sys)> BatchDerivFunc
Integrate N independent 2-equation ODE systems in batch.
Definition OdeSolver.hpp:112
constexpr int ODE_UNDERFLOW
Step size underflowed.
Definition OdeSolver.hpp:52
constexpr int ODE_OK
Success.
Definition OdeSolver.hpp:50
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.
Definition OdeSolver.cpp:114
constexpr int ODE_MAX_STEPS
Exceeded MAXSTP.
Definition OdeSolver.hpp:53
constexpr double TINY
Underflow protection.
Definition OdeSolver.hpp:40
int integrate_batch_2eq(double *y0, double *y1, int n_sys, double x1, double x2, double eps, double h1, const BatchDerivFunc &derivs)
Definition OdeSolver.cpp:146
constexpr double ERRCON
= (5/SAFETY)^(1/PGROW)
Definition OdeSolver.hpp:44
constexpr int ODE_TOO_MANY
n > allocated max
Definition OdeSolver.hpp:51
Definition Controls.cpp:24
double * dydx
Definition odesolve.c:32
double * y
Definition odesolve.c:28