Skip to main content

numra_dde/
system.rs

1//! DDE system trait and solver infrastructure.
2//!
3//! Author: Moussa Leblouba
4//! Date: 4 February 2026
5//! Modified: 2 May 2026
6
7use numra_core::Scalar;
8
9/// Trait for delay differential equation systems.
10///
11/// Defines a DDE of the form:
12/// ```text
13/// y'(t) = f(t, y(t), y(t - τ₁), y(t - τ₂), ...)
14/// ```
15pub trait DdeSystem<S: Scalar> {
16    /// Dimension of the state space.
17    fn dim(&self) -> usize;
18
19    /// Get the delays τᵢ.
20    ///
21    /// For constant delays, this returns a fixed vector.
22    /// For state-dependent delays, implement `delays_at` instead.
23    fn delays(&self) -> Vec<S>;
24
25    /// Get delays at a specific state (for state-dependent delays).
26    ///
27    /// Default implementation returns constant delays.
28    fn delays_at(&self, _t: S, _y: &[S]) -> Vec<S> {
29        self.delays()
30    }
31
32    /// Number of delays.
33    fn n_delays(&self) -> usize {
34        self.delays().len()
35    }
36
37    /// Evaluate the right-hand side f(t, y(t), y(t-τ₁), ...).
38    ///
39    /// # Arguments
40    /// * `t` - Current time
41    /// * `y` - Current state y(t)
42    /// * `y_delayed` - Delayed states: `y_delayed[i] = y(t - τᵢ)`
43    /// * `dydt` - Output buffer for y'(t)
44    fn rhs(&self, t: S, y: &[S], y_delayed: &[&[S]], dydt: &mut [S]);
45
46    /// Whether delays are state-dependent.
47    fn has_state_dependent_delays(&self) -> bool {
48        false
49    }
50}
51
52/// Options for DDE solvers.
53#[derive(Clone, Debug)]
54pub struct DdeOptions<S: Scalar> {
55    /// Relative tolerance
56    pub rtol: S,
57    /// Absolute tolerance
58    pub atol: S,
59    /// Initial step size (None = auto)
60    pub h0: Option<S>,
61    /// Maximum step size
62    pub h_max: S,
63    /// Minimum step size
64    pub h_min: S,
65    /// Maximum number of steps
66    pub max_steps: usize,
67    /// Save solution at these times (None = save all steps)
68    pub t_eval: Option<Vec<S>>,
69    /// Enable dense output (required for accurate history interpolation)
70    pub dense_output: bool,
71    /// Track discontinuities in the solution
72    pub track_discontinuities: bool,
73    /// Track discontinuities up to this derivative order
74    pub discontinuity_order: usize,
75}
76
77impl<S: Scalar> Default for DdeOptions<S> {
78    fn default() -> Self {
79        Self {
80            rtol: S::from_f64(1e-6),
81            atol: S::from_f64(1e-9),
82            h0: None,
83            h_max: S::INFINITY,
84            h_min: S::from_f64(1e-14),
85            max_steps: 100_000,
86            t_eval: None,
87            dense_output: true, // Usually needed for DDEs
88            track_discontinuities: false,
89            discontinuity_order: 0,
90        }
91    }
92}
93
94impl<S: Scalar> DdeOptions<S> {
95    pub fn rtol(mut self, rtol: S) -> Self {
96        self.rtol = rtol;
97        self
98    }
99
100    pub fn atol(mut self, atol: S) -> Self {
101        self.atol = atol;
102        self
103    }
104
105    pub fn h0(mut self, h0: S) -> Self {
106        self.h0 = Some(h0);
107        self
108    }
109
110    pub fn h_max(mut self, h_max: S) -> Self {
111        self.h_max = h_max;
112        self
113    }
114
115    pub fn max_steps(mut self, max_steps: usize) -> Self {
116        self.max_steps = max_steps;
117        self
118    }
119
120    /// Enable discontinuity tracking.
121    ///
122    /// When enabled, the solver will detect and step exactly to
123    /// discontinuity points that propagate from t0 at each delay interval.
124    pub fn track_discontinuities(mut self, track: bool) -> Self {
125        self.track_discontinuities = track;
126        self
127    }
128
129    /// Set the discontinuity order to track.
130    ///
131    /// Discontinuities propagate: the initial discontinuity at t0 propagates
132    /// to t0 + tau for each delay tau. Setting order=3 tracks up to 3 levels
133    /// of propagation (e.g., t0, t0+tau, t0+2*tau, t0+3*tau for single delay).
134    pub fn discontinuity_order(mut self, order: usize) -> Self {
135        self.discontinuity_order = order;
136        self
137    }
138}
139
140/// Statistics from DDE solver.
141#[derive(Clone, Debug, Default)]
142pub struct DdeStats {
143    /// Number of RHS evaluations
144    pub n_eval: usize,
145    /// Number of accepted steps
146    pub n_accept: usize,
147    /// Number of rejected steps
148    pub n_reject: usize,
149    /// Number of discontinuity crossings handled
150    pub n_discontinuities: usize,
151}
152
153/// Result of DDE integration.
154#[derive(Clone, Debug)]
155pub struct DdeResult<S: Scalar> {
156    /// Time points
157    pub t: Vec<S>,
158    /// Solution at each time point (row-major: y[i*dim + j] = y_j(t_i))
159    pub y: Vec<S>,
160    /// Dimension of the system
161    pub dim: usize,
162    /// Solver statistics
163    pub stats: DdeStats,
164    /// Was integration successful?
165    pub success: bool,
166    /// Message (error description if failed)
167    pub message: String,
168}
169
170impl<S: Scalar> DdeResult<S> {
171    pub fn new(t: Vec<S>, y: Vec<S>, dim: usize, stats: DdeStats) -> Self {
172        Self {
173            t,
174            y,
175            dim,
176            stats,
177            success: true,
178            message: String::new(),
179        }
180    }
181
182    pub fn failed(message: String, stats: DdeStats) -> Self {
183        Self {
184            t: Vec::new(),
185            y: Vec::new(),
186            dim: 0,
187            stats,
188            success: false,
189            message,
190        }
191    }
192
193    pub fn len(&self) -> usize {
194        self.t.len()
195    }
196
197    pub fn is_empty(&self) -> bool {
198        self.t.is_empty()
199    }
200
201    pub fn t_final(&self) -> Option<S> {
202        self.t.last().copied()
203    }
204
205    pub fn y_final(&self) -> Option<Vec<S>> {
206        if self.t.is_empty() {
207            None
208        } else {
209            let start = (self.t.len() - 1) * self.dim;
210            Some(self.y[start..start + self.dim].to_vec())
211        }
212    }
213
214    pub fn y_at(&self, i: usize) -> &[S] {
215        let start = i * self.dim;
216        &self.y[start..start + self.dim]
217    }
218
219    pub fn iter(&self) -> impl Iterator<Item = (S, &[S])> {
220        self.t
221            .iter()
222            .enumerate()
223            .map(move |(i, &t)| (t, self.y_at(i)))
224    }
225}
226
227/// Trait for DDE solvers.
228pub trait DdeSolver<S: Scalar> {
229    /// Solve the DDE problem.
230    ///
231    /// # Arguments
232    /// * `system` - The DDE system to solve
233    /// * `t0` - Initial time
234    /// * `tf` - Final time
235    /// * `history` - History function φ(t) for t ≤ t0
236    /// * `options` - Solver options
237    fn solve<Sys, H>(
238        system: &Sys,
239        t0: S,
240        tf: S,
241        history: &H,
242        options: &DdeOptions<S>,
243    ) -> Result<DdeResult<S>, String>
244    where
245        Sys: DdeSystem<S>,
246        H: Fn(S) -> Vec<S>;
247}
248
249#[cfg(test)]
250mod tests {
251    use super::*;
252
253    struct TestDde;
254
255    impl DdeSystem<f64> for TestDde {
256        fn dim(&self) -> usize {
257            1
258        }
259        fn delays(&self) -> Vec<f64> {
260            vec![1.0]
261        }
262        fn rhs(&self, _t: f64, y: &[f64], y_delayed: &[&[f64]], dydt: &mut [f64]) {
263            dydt[0] = -y[0] + y_delayed[0][0];
264        }
265    }
266
267    #[test]
268    fn test_dde_system_trait() {
269        let sys = TestDde;
270        assert_eq!(sys.dim(), 1);
271        assert_eq!(sys.n_delays(), 1);
272        assert_eq!(sys.delays(), vec![1.0]);
273
274        let mut dydt = [0.0];
275        let y_delayed = [&[0.5][..]];
276        sys.rhs(0.0, &[1.0], &y_delayed, &mut dydt);
277        assert!((dydt[0] - (-0.5)).abs() < 1e-10);
278    }
279
280    #[test]
281    fn test_dde_options() {
282        let opts: DdeOptions<f64> = DdeOptions::default().rtol(1e-8).atol(1e-10);
283        assert!((opts.rtol - 1e-8).abs() < 1e-15);
284        assert!((opts.atol - 1e-10).abs() < 1e-15);
285    }
286}