File diffeq.hpp
↰ Parent directory (include
)
Modern C++ ODE Integration Library with Real-time Signal Processing.
Definition (include/diffeq.hpp
)
Detailed Description
This library provides a comprehensive C++20 implementation of ODE integrators with advanced real-time capabilities for financial and robotics applications. Features include real-time signal processing, inter-process communication, and async execution using modern C++ standards.
Core Integrators:
Fixed Step Methods:
EulerIntegrator: Simple 1st order explicit method
ImprovedEulerIntegrator: 2nd order Heun’s method
RK4Integrator: Classic 4th order Runge-Kutta
Adaptive Step Methods (Non-stiff):
RK23Integrator: 3rd order Bogacki-Shampine with error control
RK45Integrator: 5th order Dormand-Prince (scipy’s default)
DOP853Integrator: 8th order high-accuracy method
Stiff System Methods:
RadauIntegrator: 5th order implicit Runge-Kutta
BDFIntegrator: Variable order (1-6) backward differentiation
Automatic Methods:
LSODAIntegrator: Automatic switching between Adams and BDF
SDE (Stochastic Differential Equation) Integrators:
Basic Methods:
EulerMaruyamaIntegrator: Basic SDE solver (strong order 0.5)
MilsteinIntegrator: Higher-order method with Lévy area (strong order 1.0)
SRI1Integrator: Stochastic Runge-Kutta method (strong order 1.0)
ImplicitEulerMaruyamaIntegrator: Implicit method for stiff SDEs
Advanced High-Order Methods (Strong Order 1.5):
SRAIntegrator: Stochastic Runge-Kutta for additive noise SDEs
SRIIntegrator: Stochastic Runge-Kutta for general Itô SDEs
SOSRAIntegrator: Stability-optimized SRA, robust to stiffness
SOSRIIntegrator: Stability-optimized SRI, robust to stiffness
Pre-configured Methods:
SRA1, SRA2: Different tableau variants for additive noise
SRIW1: Weak order 2.0 variant for general SDEs
SOSRA, SOSRI: Stability-optimized for high tolerances
Real-time Capabilities:
AsyncIntegrator: Lightweight async wrapper using std::future and std::thread
SignalProcessor: Generic signal processing with type-safe handlers
IntegrationInterface: Unified interface for signal-aware ODE integration
Extensible design: Works for finance, robotics, science, and any domain
Key Features:
Header-only design (no external dependencies)
Standard C++20 facilities only with proper concepts
Optional external library support (networking, JSON, etc.)
Thread-safe async operations
Unified interface for all application domains
Usage Examples:
Basic ODE Integration: #include<File diffeq.hpp> #include<vector>
//DefineODEsystem:dy/dt=-y voidexponential_decay(doublet,conststd::vector<double>&y,std::vector<double>&dydt){ dydt[0]=-y[0]; }
intmain(){ std::vector<double>y={1.0};//Initialcondition RK45Integrator<std::vector<double>>integrator(exponential_decay); integrator.integrate(y,0.1,1.0);//Integratefromt=0tot=1 //Result:y[0]≈exp(-1)≈0.368 return0; }
Real-time Signal-Aware Integration:
#include <File diffeq.hpp> #include <interfaces/integration_interface.hpp>
// Create signal-aware interface auto interface = diffeq::interfaces::make_integration_interface<std::vector<double>>();
// Register signal influences interface->register_signal_influence<double>(“price_update”, diffeq::interfaces::IntegrationInterface<std::vector<double>>::InfluenceMode::CONTINUOUS_SHIFT, [](const double& price, auto& state, auto t) { // Modify portfolio dynamics based on price signal double momentum = (price > 100.0) ? 0.01 : -0.01; for (auto& asset : state) asset *= (1.0 + momentum); });
// Register real-time output interface->register_output_stream(“monitor”, [](const auto& state, auto t) { std::cout << “Portfolio value: “ << std::accumulate(state.begin(), state.end(), 0.0) << std::endl; });
// Create signal-aware ODE auto signal_ode = interface->make_signal_aware_ode(my_portfolio_ode); auto integrator = diffeq::make_rk45<std::vector<double>>(signal_ode);
// Integration automatically handles signals and outputs integrator.integrate(state, dt, t_final);
Real-time Robot Control:
#include <File diffeq.hpp> #include <interfaces/integration_interface.hpp>
constexpr size_t N_JOINTS = 6; std::array<double, N_JOINTS * 3> robot_state{}; // position, velocity, acceleration
// Create robotics interface auto interface = diffeq::interfaces::make_integration_interface<std::array<double, 18>>();
// Register control signal influence interface->register_signal_influence<std::vector<double>>(“control_targets”, diffeq::interfaces::IntegrationInterface<std::array<double, 18>>::InfluenceMode::DISCRETE_EVENT, [](const auto& targets, auto& state, auto t) { // Update target positions for each joint for (size_t i = 0; i < targets.size() && i < N_JOINTS; ++i) { // Apply control logic } });
// Emergency stop capability interface->register_signal_influence<bool>(“emergency_stop”, diffeq::interfaces::IntegrationInterface<std::array<double, 18>>::InfluenceMode::DISCRETE_EVENT, [](bool stop, auto& state, auto t) { if (stop) { // Set all velocities to zero for (size_t i = N_JOINTS; i < 2 * N_JOINTS; ++i) state[i] = 0.0; } });
// Create robotics interface and register signals for real-time control. // Example shows emergency stop and joint monitoring capabilities.
Stochastic Differential Equations (SDEs):
#include<File diffeq.hpp> #include<sde/sde_base.hpp> #include<sde/sde_solvers.hpp>
usingnamespacediffeq::sde;
//DefineSDEsystem:dX=f(t,X)dt+g(t,X)dW //Example:GeometricBrownianMotiondS=μSdt+σSdW autodrift=[](doublet,conststd::vector<double>&x,std::vector<double>&fx){ doublemu=0.05;//5%drift fx[0]=mu*x[0]; };
autodiffusion=[](doublet,conststd::vector<double>&x,std::vector<double>&gx){ doublesigma=0.2;//20%volatility gx[0]=sigma*x[0]; };
//CreateSDEproblemandintegrator autoproblem=factory::make_sde_problem<std::vector<double>,double>( drift,diffusion,NoiseType::DIAGONAL_NOISE); autowiener=factory::make_wiener_process<std::vector<double>,double>(1,42);
Template Class EulerMaruyamaIntegrator,double>integrator(problem,wiener);
std::vector<double>state={100.0};//Initialstockprice integrator.integrate(state,0.01,1.0);//Integratetot=1
Available SDE methods:
EulerMaruyamaIntegrator: Basic method (strong order 0.5)
MilsteinIntegrator: Higher accuracy with Lévy area (strong order 1.0)
SRI1Integrator: Stochastic Runge-Kutta (strong order 1.0)
ImplicitEulerMaruyamaIntegrator: For stiff SDEs
Guidelines for Integrator Selection:
RK45Integrator: Best general-purpose choice for smooth, non-stiff ODEs. Similar to scipy.integrate.solve_ivp with method RK45.
RK23Integrator: Lower accuracy but faster for less demanding problems. Similar to scipy.integrate.solve_ivp with method RK23.
DOP853Integrator: High accuracy for demanding smooth problems. Similar to scipy.integrate.solve_ivp with method DOP853.
BDFIntegrator: Variable order method for stiff systems. Similar to scipy.integrate.solve_ivp with method BDF.
RadauIntegrator: Implicit method for stiff systems. Similar to scipy.integrate.solve_ivp with method Radau.
LSODAIntegrator: Automatic stiffness detection and method switching. Similar to scipy.integrate.odeint or solve_ivp with method LSODA.
RK4Integrator: Simple fixed-step method for educational purposes or when step size control is handled externally.
Concepts:
system_state: Container types (std::vector, std::array, etc.) that can represent the state vector of the ODE system.
can_be_time: Arithmetic types (double, float, int) that can represent time.
All integrators and interfaces are templated on these concepts for maximum flexibility and type safety.
Architecture:
The library uses a unified interface design where all real-time signal processing capabilities are provided through a single IntegrationInterface class. This interface supports:
Discrete Events: Instantaneous state modifications triggered by signals
Continuous Influences: Ongoing trajectory modifications from signals
Parameter Updates: Dynamic changes to integration parameters
Real-time Output: Streaming integration data at specified intervals
This design eliminates the need for domain-specific processors while remaining flexible enough to handle any application domain.
Includes
async/async_integrator.hpp
(File async_integrator.hpp)core/abstract_integrator.hpp
(File abstract_integrator.hpp)core/adaptive_integrator.hpp
(File adaptive_integrator.hpp)core/composable_integration.hpp
(File composable_integration.hpp)core/concepts.hpp
(File concepts.hpp)core/timeout_integrator.hpp
(File timeout_integrator.hpp)integrators/ode/bdf.hpp
(File bdf.hpp)integrators/ode/dop853.hpp
(File dop853.hpp)integrators/ode/euler.hpp
(File euler.hpp)integrators/ode/improved_euler.hpp
(File improved_euler.hpp)integrators/ode/lsoda.hpp
(File lsoda.hpp)integrators/ode/rk23.hpp
(File rk23.hpp)integrators/ode/rk4.hpp
(File rk4.hpp)integrators/ode/rk45.hpp
(File rk45.hpp)integrators/sde/euler_maruyama.hpp
(File euler_maruyama.hpp)integrators/sde/implicit_euler_maruyama.hpp
(File implicit_euler_maruyama.hpp)integrators/sde/milstein.hpp
(File milstein.hpp)integrators/sde/sosra.hpp
(File sosra.hpp)integrators/sde/sosri.hpp
(File sosri.hpp)integrators/sde/sra.hpp
(File sra.hpp)integrators/sde/sra1.hpp
(File sra1.hpp)integrators/sde/sra2.hpp
(File sra2.hpp)integrators/sde/sri.hpp
(File sri.hpp)integrators/sde/sri1.hpp
(File sri1.hpp)integrators/sde/sriw1.hpp
(File sriw1.hpp)interfaces/integration_interface.hpp
(File integration_interface.hpp)sde/sde_base.hpp
(File sde_base.hpp)signal/signal_processor.hpp
(File signal_processor.hpp)