Skip to main content

numra_spde/
solver.rs

1//! SPDE solvers using Method of Lines + SDE steppers.
2//!
3//! ## Approach
4//!
5//! Spatial discretization converts the SPDE into a system of SDEs (one per
6//! grid point), which are then integrated with an SDE stepper (Euler-Maruyama
7//! or Milstein).
8//!
9//! ## Adaptive Stepping
10//!
11//! Adaptive time stepping uses step-doubling error estimation: each step is
12//! taken once with step size h and twice with h/2, using independent Wiener
13//! increments. The difference estimates the local error.
14//!
15//! **Limitations:**
16//! - The weak convergence order of Euler-Maruyama is 1.0, so the step size
17//!   controller uses a conservative exponent of 1/3.
18//! - Adaptive stepping adds overhead from the extra half-steps and noise
19//!   generation.
20//! - For problems with very stiff diffusion operators, fixed-step methods
21//!   may be more efficient.
22//!
23//! ## Noise Correlation
24//!
25//! Three noise types are supported via [`NoiseCorrelation`]:
26//! - **White:** Independent noise at each spatial point (diagonal diffusion)
27//! - **Spatial:** Correlated noise with a given correlation length (uses
28//!   Cholesky factorization of the spatial correlation matrix)
29//! - **Trace-class:** Correlated noise with eigenvalue decay (Q-Wiener process)
30//!
31//! Author: Moussa Leblouba
32//! Date: 8 February 2026
33//! Modified: 2 May 2026
34
35use numra_core::Scalar;
36use numra_pde::Grid1D;
37use numra_sde::{
38    create_wiener, EulerMaruyama, Milstein, NoiseType, SdeOptions, SdeSolver, SdeSystem,
39};
40
41use crate::noise::cholesky_decompose;
42use crate::system::{NoiseCorrelation, SpdeSystem};
43
44/// SDE integration method for the spatially-discretized system.
45#[derive(Clone, Copy, Debug, Default)]
46pub enum SpdeMethod {
47    /// Euler-Maruyama (1st order, simple)
48    #[default]
49    EulerMaruyama,
50    /// Milstein (1.0 order strong convergence, for multiplicative noise)
51    Milstein,
52}
53
54/// Options for SPDE solver.
55///
56/// **Divergence from `numra_ode::SolverOptions`** (per Foundation Spec §2.5):
57/// SPDE solvers carry both stochastic noise (like SDE — Wiener time
58/// control, `seed`) and an explicit method selector (`SpdeMethod`,
59/// `EulerMaruyama` | `Milstein`). The method selector is unusual:
60/// `numra-ode` solvers dispatch at the type level
61/// (`DoPri5::solve(...)`, `Tsit5::solve(...)`); SPDE collapses to a single
62/// `solve` entry that branches on the enum because noise-discretisation
63/// logic is shared across methods. `n_output` discretises sampled-output
64/// count rather than `SolverOptions::t_eval`-style explicit time grid —
65/// appropriate for stochastic workloads where exact output times matter
66/// less than density. `adaptive: bool` gates between fixed-step EM
67/// (default) and adaptive variants. See
68/// `docs/architecture/foundation-specification.md` §2.5.
69#[derive(Clone, Debug)]
70pub struct SpdeOptions<S: Scalar> {
71    /// Time step
72    pub dt: S,
73    /// Number of output time points
74    pub n_output: usize,
75    /// SDE method
76    pub method: SpdeMethod,
77    /// Random seed
78    pub seed: Option<u64>,
79    /// Enable adaptive time stepping
80    pub adaptive: bool,
81    /// Relative tolerance for adaptive stepping
82    pub rtol: S,
83    /// Absolute tolerance for adaptive stepping
84    pub atol: S,
85    /// Minimum time step for adaptive stepping
86    pub dt_min: S,
87    /// Maximum time step for adaptive stepping
88    pub dt_max: S,
89}
90
91impl<S: Scalar> Default for SpdeOptions<S> {
92    fn default() -> Self {
93        Self {
94            dt: S::from_f64(0.001),
95            n_output: 100,
96            method: SpdeMethod::EulerMaruyama,
97            seed: None,
98            adaptive: false,
99            rtol: S::from_f64(1e-3),
100            atol: S::from_f64(1e-6),
101            dt_min: S::from_f64(1e-10),
102            dt_max: S::from_f64(0.1),
103        }
104    }
105}
106
107impl<S: Scalar> SpdeOptions<S> {
108    /// Set time step.
109    pub fn dt(mut self, dt: S) -> Self {
110        self.dt = dt;
111        self
112    }
113
114    /// Set number of output points.
115    pub fn n_output(mut self, n: usize) -> Self {
116        self.n_output = n;
117        self
118    }
119
120    /// Set SDE method.
121    pub fn method(mut self, method: SpdeMethod) -> Self {
122        self.method = method;
123        self
124    }
125
126    /// Set random seed.
127    pub fn seed(mut self, seed: u64) -> Self {
128        self.seed = Some(seed);
129        self
130    }
131
132    /// Enable/disable adaptive time stepping.
133    pub fn adaptive(mut self, enable: bool) -> Self {
134        self.adaptive = enable;
135        self
136    }
137
138    /// Set relative tolerance for adaptive stepping.
139    pub fn rtol(mut self, tol: S) -> Self {
140        self.rtol = tol;
141        self
142    }
143
144    /// Set absolute tolerance for adaptive stepping.
145    pub fn atol(mut self, tol: S) -> Self {
146        self.atol = tol;
147        self
148    }
149
150    /// Set minimum time step for adaptive stepping.
151    pub fn dt_min(mut self, dt: S) -> Self {
152        self.dt_min = dt;
153        self
154    }
155
156    /// Set maximum time step for adaptive stepping.
157    pub fn dt_max(mut self, dt: S) -> Self {
158        self.dt_max = dt;
159        self
160    }
161}
162
163/// Statistics for SPDE solve.
164#[derive(Clone, Debug, Default)]
165pub struct SpdeStats {
166    /// Number of time steps taken
167    pub n_steps: usize,
168    /// Number of drift evaluations
169    pub n_drift: usize,
170    /// Number of diffusion evaluations
171    pub n_diffusion: usize,
172    /// Number of rejected steps (for adaptive methods)
173    pub n_reject: usize,
174}
175
176/// Result of SPDE solve.
177#[derive(Clone, Debug)]
178pub struct SpdeResult<S: Scalar> {
179    /// Time points
180    pub t: Vec<S>,
181    /// Solution at each time: `y[i]` is solution at `t[i]`
182    pub y: Vec<Vec<S>>,
183    /// Spatial grid
184    pub grid: Grid1D<S>,
185    /// Solver statistics
186    pub stats: SpdeStats,
187    /// Whether solve was successful
188    pub success: bool,
189}
190
191impl<S: Scalar> SpdeResult<S> {
192    /// Get solution at final time.
193    pub fn y_final(&self) -> Option<&[S]> {
194        self.y.last().map(|v| v.as_slice())
195    }
196
197    /// Get final time.
198    pub fn t_final(&self) -> Option<S> {
199        self.t.last().copied()
200    }
201
202    /// Get solution at time index i.
203    pub fn y_at(&self, i: usize) -> &[S] {
204        &self.y[i]
205    }
206}
207
208/// Trait for SPDE solvers.
209pub trait SpdeSolver<S: Scalar> {
210    /// Solve the SPDE.
211    fn solve<Sys: SpdeSystem<S> + Sync>(
212        system: &Sys,
213        t0: S,
214        tf: S,
215        u0: &[S],
216        grid: &Grid1D<S>,
217        options: &SpdeOptions<S>,
218    ) -> Result<SpdeResult<S>, String>;
219}
220
221/// Pre-computed noise configuration for the MOL wrapper.
222enum MolNoiseConfig<S: Scalar> {
223    /// White (diagonal) noise — each grid point gets independent noise.
224    White,
225    /// Colored noise with pre-computed Cholesky factor L (dim × dim, row-major).
226    /// Correlated noise = L * independent_noise.
227    Colored { cholesky: Vec<S>, dim: usize },
228    /// Trace-class noise using spectral representation.
229    /// n_modes independent Wiener processes drive the spatial modes.
230    TraceClass {
231        /// sqrt(lambda_k) for each mode k
232        sqrt_eigenvalues: Vec<S>,
233        /// Basis function values: basis[i * n_modes + k] = e_k(x_i)
234        basis: Vec<S>,
235        n_modes: usize,
236        dim: usize,
237    },
238}
239
240impl<S: Scalar> MolNoiseConfig<S> {
241    /// Build noise configuration from SPDE system and grid.
242    fn from_system<Sys: SpdeSystem<S>>(system: &Sys, grid: &Grid1D<S>) -> Result<Self, String> {
243        match system.noise_correlation() {
244            NoiseCorrelation::White => Ok(MolNoiseConfig::White),
245            NoiseCorrelation::Colored { correlation_length } => {
246                let dim = grid.len() - 2; // interior points
247                let dx = grid.dx_uniform();
248
249                // Build correlation matrix: C_ij = exp(-|x_i - x_j| / lambda)
250                let mut c = vec![S::ZERO; dim * dim];
251                for i in 0..dim {
252                    for j in 0..dim {
253                        let dist = dx * S::from_usize(if i > j { i - j } else { j - i });
254                        c[i * dim + j] = (-dist / correlation_length).exp();
255                    }
256                }
257
258                let cholesky = cholesky_decompose(&c, dim)?;
259                Ok(MolNoiseConfig::Colored { cholesky, dim })
260            }
261            NoiseCorrelation::TraceClass {
262                n_modes,
263                decay_rate,
264            } => {
265                let dim = grid.len() - 2;
266                let interior = grid.interior_points();
267                let (x_min, x_max) = grid.bounds();
268                let domain_length = x_max - x_min;
269
270                // Eigenvalues: lambda_k = k^(-decay_rate)
271                let sqrt_eigenvalues: Vec<S> = (1..=n_modes)
272                    .map(|k| S::from_usize(k).powf(-decay_rate).sqrt())
273                    .collect();
274
275                // Basis functions: e_k(x) = sqrt(2/L) * sin(k*pi*x/L)
276                let norm = (S::TWO / domain_length).sqrt();
277                let pi = S::PI;
278                let mut basis = vec![S::ZERO; dim * n_modes];
279                for i in 0..dim {
280                    let x = interior[i] - x_min;
281                    for k in 0..n_modes {
282                        let mode = S::from_usize(k + 1);
283                        basis[i * n_modes + k] = norm * (mode * pi * x / domain_length).sin();
284                    }
285                }
286
287                Ok(MolNoiseConfig::TraceClass {
288                    sqrt_eigenvalues,
289                    basis,
290                    n_modes,
291                    dim,
292                })
293            }
294        }
295    }
296
297    /// Number of Wiener processes needed.
298    fn n_wiener(&self, dim: usize) -> usize {
299        match self {
300            MolNoiseConfig::White => dim,
301            MolNoiseConfig::Colored { dim: d, .. } => *d,
302            MolNoiseConfig::TraceClass { n_modes, .. } => *n_modes,
303        }
304    }
305}
306
307/// Wrapper that converts SPDE to SDE system via Method of Lines.
308struct MolSdeWrapper<'a, S: Scalar, Sys: SpdeSystem<S> + Sync> {
309    spde: &'a Sys,
310    grid: &'a Grid1D<S>,
311    noise_config: MolNoiseConfig<S>,
312}
313
314impl<'a, S: Scalar, Sys: SpdeSystem<S> + Sync> SdeSystem<S> for MolSdeWrapper<'a, S, Sys> {
315    fn dim(&self) -> usize {
316        self.grid.len() - 2 // Interior points only
317    }
318
319    fn noise_type(&self) -> NoiseType {
320        match &self.noise_config {
321            MolNoiseConfig::White => NoiseType::Diagonal,
322            MolNoiseConfig::Colored { dim, .. } => NoiseType::General { n_wiener: *dim },
323            MolNoiseConfig::TraceClass { n_modes, .. } => NoiseType::General { n_wiener: *n_modes },
324        }
325    }
326
327    fn drift(&self, t: S, y: &[S], f: &mut [S]) {
328        self.spde.drift(t, y, f, self.grid);
329    }
330
331    fn diffusion(&self, t: S, y: &[S], g: &mut [S]) {
332        match &self.noise_config {
333            MolNoiseConfig::White => {
334                // Diagonal: g[i] = sigma(u_i)
335                self.spde.diffusion(t, y, g, self.grid);
336            }
337            MolNoiseConfig::Colored { cholesky, dim } => {
338                // General: g[i * dim + j] = sigma(u_i) * L[i][j]
339                // Get the base diffusion coefficients into a local buffer
340                let mut sigma = vec![S::ZERO; *dim];
341                self.spde.diffusion(t, y, &mut sigma, self.grid);
342
343                // Build g[i * dim + j] = sigma[i] * L[i][j]
344                for i in 0..*dim {
345                    for j in 0..=i {
346                        g[i * dim + j] = sigma[i] * cholesky[i * dim + j];
347                    }
348                    // Upper triangle is zero (L is lower triangular)
349                    for j in (i + 1)..*dim {
350                        g[i * dim + j] = S::ZERO;
351                    }
352                }
353            }
354            MolNoiseConfig::TraceClass {
355                sqrt_eigenvalues,
356                basis,
357                n_modes,
358                dim,
359            } => {
360                // General: g[i * n_modes + k] = sigma(u_i) * sqrt(lambda_k) * e_k(x_i)
361                let mut sigma = vec![S::ZERO; *dim];
362                self.spde.diffusion(t, y, &mut sigma, self.grid);
363
364                for i in 0..*dim {
365                    for k in 0..*n_modes {
366                        g[i * n_modes + k] =
367                            sigma[i] * sqrt_eigenvalues[k] * basis[i * n_modes + k];
368                    }
369                }
370            }
371        }
372    }
373}
374
375// SAFETY: MolSdeWrapper contains only shared references (&'a Sys where Sys: Sync,
376// &'a Grid1D<S>) and owned data (MolNoiseConfig<S> which is Send+Sync since S: Scalar
377// is Send+Sync and Vec<S> is Send+Sync). The Sys: Sync bound ensures the shared system
378// reference is safe to access from multiple threads. MolNoiseConfig only contains
379// Vec<S> and usize fields, which are inherently Sync. The compiler cannot auto-derive
380// Sync because of the generic Sys parameter interaction with the lifetime, but all
381// fields are safe to share across threads.
382unsafe impl<'a, S: Scalar, Sys: SpdeSystem<S> + Sync> Sync for MolSdeWrapper<'a, S, Sys> {}
383
384/// Method of Lines SPDE solver.
385///
386/// Converts the SPDE to a large SDE system by discretizing spatial derivatives,
387/// then applies an SDE solver.
388pub struct MolSdeSolver;
389
390impl<S: Scalar> SpdeSolver<S> for MolSdeSolver {
391    fn solve<Sys: SpdeSystem<S> + Sync>(
392        system: &Sys,
393        t0: S,
394        tf: S,
395        u0: &[S],
396        grid: &Grid1D<S>,
397        options: &SpdeOptions<S>,
398    ) -> Result<SpdeResult<S>, String> {
399        let n_interior = grid.len() - 2;
400        if u0.len() != n_interior {
401            return Err(format!(
402                "Initial condition length {} doesn't match interior points {}",
403                u0.len(),
404                n_interior
405            ));
406        }
407
408        // Use adaptive or fixed time stepping
409        if options.adaptive {
410            Self::solve_adaptive(system, t0, tf, u0, grid, options)
411        } else {
412            Self::solve_fixed(system, t0, tf, u0, grid, options)
413        }
414    }
415}
416
417impl MolSdeSolver {
418    /// Solve using fixed time stepping (delegates to SDE solver).
419    fn solve_fixed<S: Scalar, Sys: SpdeSystem<S> + Sync>(
420        system: &Sys,
421        t0: S,
422        tf: S,
423        u0: &[S],
424        grid: &Grid1D<S>,
425        options: &SpdeOptions<S>,
426    ) -> Result<SpdeResult<S>, String> {
427        // Create MOL SDE wrapper with pre-computed noise configuration
428        let noise_config = MolNoiseConfig::from_system(system, grid)?;
429        let sde_system = MolSdeWrapper {
430            spde: system,
431            grid,
432            noise_config,
433        };
434
435        // SDE options
436        let sde_options = SdeOptions::default().dt(options.dt);
437
438        // Solve using appropriate method
439        let sde_result = match options.method {
440            SpdeMethod::EulerMaruyama => {
441                EulerMaruyama::solve(&sde_system, t0, tf, u0, &sde_options, options.seed)
442            }
443            SpdeMethod::Milstein => {
444                if system.is_additive() {
445                    // For additive noise, Milstein reduces to Euler-Maruyama
446                    EulerMaruyama::solve(&sde_system, t0, tf, u0, &sde_options, options.seed)
447                } else {
448                    Milstein::solve(&sde_system, t0, tf, u0, &sde_options, options.seed)
449                }
450            }
451        }?;
452
453        // Build output with subsampling
454        let n_steps = sde_result.t.len();
455        let stride = (n_steps / options.n_output).max(1);
456
457        let mut t_out = Vec::new();
458        let mut y_out = Vec::new();
459
460        for i in (0..n_steps).step_by(stride) {
461            t_out.push(sde_result.t[i]);
462            y_out.push(sde_result.y_at(i).to_vec());
463        }
464
465        // Always include final point
466        if t_out.last() != sde_result.t.last() {
467            t_out.push(*sde_result.t.last().unwrap());
468            y_out.push(sde_result.y_final().unwrap().to_vec());
469        }
470
471        Ok(SpdeResult {
472            t: t_out,
473            y: y_out,
474            grid: grid.clone(),
475            stats: SpdeStats {
476                n_steps: sde_result.stats.n_accept + sde_result.stats.n_reject,
477                n_drift: sde_result.stats.n_drift,
478                n_diffusion: sde_result.stats.n_diffusion,
479                n_reject: sde_result.stats.n_reject,
480            },
481            success: sde_result.success,
482        })
483    }
484
485    /// Solve using adaptive time stepping with step doubling for error estimation.
486    fn solve_adaptive<S: Scalar, Sys: SpdeSystem<S> + Sync>(
487        system: &Sys,
488        t0: S,
489        tf: S,
490        u0: &[S],
491        grid: &Grid1D<S>,
492        options: &SpdeOptions<S>,
493    ) -> Result<SpdeResult<S>, String> {
494        let dim = u0.len();
495        let mut t = t0;
496        let mut u = u0.to_vec();
497        let mut dt = options.dt;
498
499        // Pre-compute noise configuration for correlated noise
500        let noise_config = MolNoiseConfig::from_system(system, grid)?;
501        let n_wiener = noise_config.n_wiener(dim);
502
503        // Wiener process for generating noise (independent increments)
504        let mut wiener = create_wiener(n_wiener, options.seed);
505
506        // Storage for results
507        let mut t_out = vec![t0];
508        let mut y_out = vec![u0.to_vec()];
509        let mut stats = SpdeStats::default();
510
511        // Buffers for computation
512        let mut drift = vec![S::ZERO; dim];
513        let mut diffusion = vec![S::ZERO; dim];
514        let mut u1 = vec![S::ZERO; dim]; // One full step
515        let mut u_half = vec![S::ZERO; dim]; // After first half step
516        let mut u2 = vec![S::ZERO; dim]; // Two half steps
517
518        let max_steps = 1_000_000;
519        let mut step_count = 0;
520
521        while t < tf && step_count < max_steps {
522            // Clamp dt to not overshoot final time
523            let h = dt.min(tf - t).min(options.dt_max).max(options.dt_min);
524
525            // Generate two independent half-step Wiener increments.
526            // dW1 ~ N(0, h/2), dW2 ~ N(0, h/2), so dW1 + dW2 ~ N(0, h).
527            let half_h = h / S::from_f64(2.0);
528            let raw_dw1 = wiener.increment(half_h);
529            let raw_dw2 = wiener.increment(half_h);
530
531            // Apply noise correlation to get effective noise at each grid point.
532            // For white noise, this is just the raw Wiener increments.
533            // For colored/trace-class, we apply the correlation structure.
534            let noise1 = apply_noise_correlation(&noise_config, &raw_dw1.dw, dim);
535            let noise2 = apply_noise_correlation(&noise_config, &raw_dw2.dw, dim);
536
537            // Compute drift and diffusion at current state
538            let mut drift_orig = vec![S::ZERO; dim];
539            let mut diffusion_orig = vec![S::ZERO; dim];
540            system.drift(t, &u, &mut drift_orig, grid);
541            system.diffusion(t, &u, &mut diffusion_orig, grid);
542            stats.n_drift += 1;
543            stats.n_diffusion += 1;
544
545            // Full step uses effective noise = noise1 + noise2
546            for i in 0..dim {
547                u1[i] = u[i] + drift_orig[i] * h + diffusion_orig[i] * (noise1[i] + noise2[i]);
548            }
549
550            // First half step uses noise1
551            for i in 0..dim {
552                u_half[i] = u[i] + drift_orig[i] * half_h + diffusion_orig[i] * noise1[i];
553            }
554
555            // Compute drift and diffusion at half step
556            system.drift(t + half_h, &u_half, &mut drift, grid);
557            system.diffusion(t + half_h, &u_half, &mut diffusion, grid);
558            stats.n_drift += 1;
559            stats.n_diffusion += 1;
560
561            // Second half step uses noise2
562            for i in 0..dim {
563                u2[i] = u_half[i] + drift[i] * half_h + diffusion[i] * noise2[i];
564            }
565
566            // Estimate error as difference between methods
567            let mut err_max = S::ZERO;
568            let mut scale_max = S::ZERO;
569            for i in 0..dim {
570                let err_i = (u1[i] - u2[i]).abs();
571                let scale_i = options.atol + options.rtol * u[i].abs().max(u2[i].abs());
572                let ratio = err_i / scale_i;
573                if ratio > err_max / scale_max.max(S::from_f64(1e-15)) {
574                    err_max = err_i;
575                    scale_max = scale_i;
576                }
577            }
578
579            let err_ratio = err_max / scale_max.max(S::from_f64(1e-15));
580
581            if err_ratio <= S::ONE || h <= options.dt_min {
582                // Accept step - use the more accurate two-step solution
583                t = t + h;
584                u.copy_from_slice(&u2);
585                stats.n_steps += 1;
586
587                // Store result
588                t_out.push(t);
589                y_out.push(u.clone());
590
591                // Increase step size if error is small
592                // Use exponent 1/3 (conservative for EM weak order 1 with step-doubling)
593                if err_ratio < S::from_f64(0.5) {
594                    let factor = (S::ONE / err_ratio.max(S::from_f64(1e-10)))
595                        .powf(S::from_f64(1.0 / 3.0))
596                        .min(S::from_f64(2.0));
597                    dt = (h * factor).min(options.dt_max);
598                }
599            } else {
600                // Reject step and decrease dt
601                stats.n_reject += 1;
602                let factor = (S::ONE / err_ratio.max(S::from_f64(1e-10)))
603                    .powf(S::from_f64(1.0 / 3.0))
604                    .max(S::from_f64(0.25))
605                    .min(S::from_f64(0.9));
606                dt = (h * factor).max(options.dt_min);
607            }
608
609            step_count += 1;
610        }
611
612        if step_count >= max_steps && t < tf {
613            return Err(format!(
614                "Maximum steps ({}) exceeded at t = {}",
615                max_steps,
616                t.to_f64()
617            ));
618        }
619
620        // Subsample output if needed
621        let n_total = t_out.len();
622        if n_total > options.n_output {
623            let stride = (n_total / options.n_output).max(1);
624            let mut t_sub = Vec::new();
625            let mut y_sub = Vec::new();
626
627            for i in (0..n_total).step_by(stride) {
628                t_sub.push(t_out[i]);
629                y_sub.push(y_out[i].clone());
630            }
631
632            // Always include final point
633            if t_sub.last() != t_out.last() {
634                t_sub.push(*t_out.last().unwrap());
635                y_sub.push(y_out.last().unwrap().clone());
636            }
637
638            t_out = t_sub;
639            y_out = y_sub;
640        }
641
642        Ok(SpdeResult {
643            t: t_out,
644            y: y_out,
645            grid: grid.clone(),
646            stats,
647            success: true,
648        })
649    }
650}
651
652/// Apply noise correlation structure to raw independent Wiener increments.
653///
654/// For white noise, returns the raw increments unchanged (identity mapping).
655/// For colored noise, applies Cholesky factor: effective_dW = L * raw_dZ.
656/// For trace-class noise, applies spectral expansion: effective_dW[i] = sum_k sqrt(lam_k) * e_k(x_i) * raw_dZ[k].
657fn apply_noise_correlation<S: Scalar>(
658    config: &MolNoiseConfig<S>,
659    raw_dw: &[S],
660    dim: usize,
661) -> Vec<S> {
662    match config {
663        MolNoiseConfig::White => {
664            // No correlation — raw increments are already per grid point
665            raw_dw[..dim].to_vec()
666        }
667        MolNoiseConfig::Colored { cholesky, dim: n } => {
668            // effective_dW[i] = sum_{j<=i} L[i][j] * raw_dZ[j]
669            let mut result = vec![S::ZERO; *n];
670            for i in 0..*n {
671                for j in 0..=i {
672                    result[i] = result[i] + cholesky[i * n + j] * raw_dw[j];
673                }
674            }
675            result
676        }
677        MolNoiseConfig::TraceClass {
678            sqrt_eigenvalues,
679            basis,
680            n_modes,
681            dim: n,
682        } => {
683            // effective_dW[i] = sum_k sqrt(lambda_k) * e_k(x_i) * raw_dZ[k]
684            let mut result = vec![S::ZERO; *n];
685            for i in 0..*n {
686                for k in 0..*n_modes {
687                    result[i] =
688                        result[i] + sqrt_eigenvalues[k] * basis[i * n_modes + k] * raw_dw[k];
689                }
690            }
691            result
692        }
693    }
694}
695
696/// Ensemble solver for SPDE Monte Carlo simulations.
697pub struct SpdeEnsemble;
698
699impl SpdeEnsemble {
700    /// Run ensemble of SPDE solutions with different noise realizations.
701    pub fn solve<S: Scalar, Sys: SpdeSystem<S> + Sync>(
702        system: &Sys,
703        t0: S,
704        tf: S,
705        u0: &[S],
706        grid: &Grid1D<S>,
707        options: &SpdeOptions<S>,
708        n_trajectories: usize,
709    ) -> Result<Vec<SpdeResult<S>>, String> {
710        let mut results = Vec::with_capacity(n_trajectories);
711        let base_seed = options.seed.unwrap_or(42);
712
713        for i in 0..n_trajectories {
714            let mut opts = options.clone();
715            opts.seed = Some(base_seed + i as u64);
716            let result = MolSdeSolver::solve(system, t0, tf, u0, grid, &opts)?;
717            results.push(result);
718        }
719
720        Ok(results)
721    }
722
723    /// Compute mean trajectory from ensemble.
724    pub fn mean<S: Scalar>(results: &[SpdeResult<S>]) -> Vec<Vec<S>> {
725        if results.is_empty() {
726            return Vec::new();
727        }
728
729        let n_times = results[0].t.len();
730        let n_points = results[0].y[0].len();
731        let n_traj = results.len();
732
733        let mut mean = vec![vec![S::ZERO; n_points]; n_times];
734
735        for result in results {
736            for (i, y) in result.y.iter().enumerate() {
737                for (j, &val) in y.iter().enumerate() {
738                    mean[i][j] = mean[i][j] + val;
739                }
740            }
741        }
742
743        let n_traj_s = S::from_usize(n_traj);
744        for m in mean.iter_mut() {
745            for val in m.iter_mut() {
746                *val = *val / n_traj_s;
747            }
748        }
749
750        mean
751    }
752
753    /// Compute standard deviation from ensemble.
754    pub fn std<S: Scalar>(results: &[SpdeResult<S>], mean: &[Vec<S>]) -> Vec<Vec<S>> {
755        if results.is_empty() {
756            return Vec::new();
757        }
758
759        let n_times = results[0].t.len();
760        let n_points = results[0].y[0].len();
761        let n_traj = results.len();
762
763        let mut variance = vec![vec![S::ZERO; n_points]; n_times];
764
765        for result in results {
766            for (i, y) in result.y.iter().enumerate() {
767                for (j, &val) in y.iter().enumerate() {
768                    let diff = val - mean[i][j];
769                    variance[i][j] = variance[i][j] + diff * diff;
770                }
771            }
772        }
773
774        let n_traj_s = S::from_usize(n_traj - 1); // Bessel's correction
775        for v in variance.iter_mut() {
776            for val in v.iter_mut() {
777                *val = (*val / n_traj_s).sqrt();
778            }
779        }
780
781        variance
782    }
783}
784
785#[cfg(test)]
786mod tests {
787    use super::*;
788
789    struct StochasticHeat {
790        alpha: f64,
791        sigma: f64,
792    }
793
794    impl SpdeSystem<f64> for StochasticHeat {
795        fn drift(&self, _t: f64, u: &[f64], du: &mut [f64], grid: &Grid1D<f64>) {
796            let dx = grid.dx_uniform();
797            let n = u.len();
798
799            for i in 0..n {
800                let u_left = if i == 0 { 0.0 } else { u[i - 1] };
801                let u_right = if i == n - 1 { 0.0 } else { u[i + 1] };
802                du[i] = self.alpha * (u_left - 2.0 * u[i] + u_right) / (dx * dx);
803            }
804        }
805
806        fn diffusion(&self, _t: f64, _u: &[f64], sigma: &mut [f64], _grid: &Grid1D<f64>) {
807            for s in sigma.iter_mut() {
808                *s = self.sigma;
809            }
810        }
811    }
812
813    #[test]
814    fn test_stochastic_heat_basic() {
815        let system = StochasticHeat {
816            alpha: 0.01,
817            sigma: 0.1,
818        };
819
820        let grid = Grid1D::uniform(0.0, 1.0, 21);
821        let u0: Vec<f64> = grid
822            .interior_points()
823            .iter()
824            .map(|&x| (std::f64::consts::PI * x).sin())
825            .collect();
826
827        let options = SpdeOptions::default().dt(0.0001).n_output(50).seed(12345);
828
829        let result = MolSdeSolver::solve(&system, 0.0, 0.1, &u0, &grid, &options).unwrap();
830
831        assert!(result.success);
832        assert!(!result.t.is_empty());
833        assert!(!result.y.is_empty());
834
835        // Solution should be finite
836        for y in &result.y {
837            for &val in y {
838                assert!(val.is_finite(), "Solution should be finite");
839            }
840        }
841    }
842
843    #[test]
844    fn test_spde_reproducibility() {
845        let system = StochasticHeat {
846            alpha: 0.01,
847            sigma: 0.1,
848        };
849
850        let grid = Grid1D::uniform(0.0, 1.0, 11);
851        let u0: Vec<f64> = grid
852            .interior_points()
853            .iter()
854            .map(|&x| (std::f64::consts::PI * x).sin())
855            .collect();
856
857        let options = SpdeOptions::default().dt(0.001).n_output(10).seed(42);
858
859        let result1 = MolSdeSolver::solve(&system, 0.0, 0.01, &u0, &grid, &options).unwrap();
860        let result2 = MolSdeSolver::solve(&system, 0.0, 0.01, &u0, &grid, &options).unwrap();
861
862        let y1 = result1.y_final().unwrap();
863        let y2 = result2.y_final().unwrap();
864
865        for (a, b) in y1.iter().zip(y2.iter()) {
866            assert!(
867                (a - b).abs() < 1e-10,
868                "Results should be reproducible with same seed"
869            );
870        }
871    }
872
873    #[test]
874    fn test_ensemble_statistics() {
875        let system = StochasticHeat {
876            alpha: 0.01,
877            sigma: 0.05,
878        };
879
880        let grid = Grid1D::uniform(0.0, 1.0, 11);
881        let u0: Vec<f64> = grid
882            .interior_points()
883            .iter()
884            .map(|&x| (std::f64::consts::PI * x).sin())
885            .collect();
886
887        let options = SpdeOptions::default().dt(0.001).n_output(10).seed(100);
888
889        let results = SpdeEnsemble::solve(&system, 0.0, 0.01, &u0, &grid, &options, 10).unwrap();
890
891        assert_eq!(results.len(), 10);
892
893        let mean = SpdeEnsemble::mean(&results);
894        let std = SpdeEnsemble::std(&results, &mean);
895
896        // Mean should be close to deterministic solution
897        // Standard deviation should be positive
898        for s in &std[std.len() - 1] {
899            assert!(*s >= 0.0, "Std should be non-negative");
900        }
901    }
902
903    #[test]
904    fn test_zero_noise() {
905        // With zero noise, should reduce to deterministic PDE
906        let system = StochasticHeat {
907            alpha: 0.1,
908            sigma: 0.0,
909        };
910
911        let grid = Grid1D::uniform(0.0, 1.0, 11);
912        let u0: Vec<f64> = grid
913            .interior_points()
914            .iter()
915            .map(|&x| (std::f64::consts::PI * x).sin())
916            .collect();
917
918        let options = SpdeOptions::default().dt(0.0001).n_output(10).seed(42);
919
920        let result1 = MolSdeSolver::solve(&system, 0.0, 0.01, &u0, &grid, &options).unwrap();
921
922        // Change seed - result should still be the same with zero noise
923        let options2 = options.seed(999);
924        let result2 = MolSdeSolver::solve(&system, 0.0, 0.01, &u0, &grid, &options2).unwrap();
925
926        let y1 = result1.y_final().unwrap();
927        let y2 = result2.y_final().unwrap();
928
929        for (a, b) in y1.iter().zip(y2.iter()) {
930            assert!(
931                (a - b).abs() < 1e-8,
932                "Zero noise should give deterministic results"
933            );
934        }
935    }
936
937    #[test]
938    fn test_options_builder() {
939        let opts: SpdeOptions<f64> = SpdeOptions::default()
940            .dt(0.01)
941            .n_output(200)
942            .method(SpdeMethod::Milstein)
943            .seed(123);
944
945        assert!((opts.dt - 0.01).abs() < 1e-10);
946        assert_eq!(opts.n_output, 200);
947        assert!(matches!(opts.method, SpdeMethod::Milstein));
948        assert_eq!(opts.seed, Some(123));
949    }
950
951    #[test]
952    fn test_spde_adaptive_dt() {
953        // Stochastic heat equation with adaptive time stepping
954        let spde = StochasticHeat {
955            alpha: 0.1,
956            sigma: 0.01,
957        };
958
959        let grid = Grid1D::uniform(0.0, 1.0, 52); // 52 total points -> 50 interior
960        let u0: Vec<f64> = grid
961            .interior_points()
962            .iter()
963            .map(|&x| (std::f64::consts::PI * x).sin())
964            .collect();
965
966        let options = SpdeOptions::default()
967            .adaptive(true)
968            .dt(0.001)
969            .rtol(1e-4)
970            .atol(1e-6)
971            .dt_max(0.01)
972            .seed(42);
973
974        let result = MolSdeSolver::solve(&spde, 0.0, 0.5, &u0, &grid, &options).unwrap();
975
976        assert!(result.success, "Adaptive solve should succeed");
977        assert!(!result.t.is_empty(), "Should have output time points");
978
979        // Solution should be finite
980        for y in &result.y {
981            for &val in y {
982                assert!(val.is_finite(), "Solution should be finite, got {}", val);
983            }
984        }
985
986        // Verify dt varies (not all equal)
987        let dts: Vec<_> = result.t.windows(2).map(|w| w[1] - w[0]).collect();
988        let first_dt = dts.first().copied().unwrap_or(0.0);
989        let has_variation = dts.iter().any(|&dt| (dt - first_dt).abs() > 1e-10);
990
991        assert!(
992            has_variation,
993            "Time steps should vary with adaptive stepping"
994        );
995    }
996
997    #[test]
998    fn test_adaptive_noise_consistency() {
999        // Verify that the adaptive solver's noise splitting is correct.
1000        // Two independent half-step increments dW1 ~ N(0, h/2), dW2 ~ N(0, h/2)
1001        // should produce a total dW1 + dW2 ~ N(0, h).
1002        //
1003        // We test this by running the adaptive solver on a pure-noise problem
1004        // (zero drift) and checking the variance of the ensemble matches theory.
1005        struct PureNoise;
1006
1007        impl SpdeSystem<f64> for PureNoise {
1008            fn drift(&self, _t: f64, _u: &[f64], du: &mut [f64], _grid: &Grid1D<f64>) {
1009                for d in du.iter_mut() {
1010                    *d = 0.0;
1011                }
1012            }
1013
1014            fn diffusion(&self, _t: f64, _u: &[f64], sigma: &mut [f64], _grid: &Grid1D<f64>) {
1015                for s in sigma.iter_mut() {
1016                    *s = 1.0;
1017                }
1018            }
1019        }
1020
1021        let n_trajectories = 500;
1022        let grid = Grid1D::uniform(0.0, 1.0, 12); // 10 interior points
1023        let n_interior = grid.len() - 2;
1024        let u0 = vec![0.0; n_interior];
1025        let tf = 0.1;
1026
1027        let mut sum_sq = vec![0.0; n_interior];
1028        for i in 0..n_trajectories {
1029            let options = SpdeOptions::default()
1030                .adaptive(true)
1031                .dt(0.01)
1032                .rtol(1e-2)
1033                .atol(1e-4)
1034                .dt_min(1e-6)
1035                .dt_max(0.05)
1036                .seed(1000 + i as u64);
1037
1038            let result = MolSdeSolver::solve(&PureNoise, 0.0, tf, &u0, &grid, &options).unwrap();
1039            let y_final = result.y_final().unwrap();
1040            for j in 0..n_interior {
1041                sum_sq[j] += y_final[j] * y_final[j];
1042            }
1043        }
1044
1045        // For pure noise with sigma=1, Var[u(T)] = T for each component
1046        // (Wiener process variance). With adaptive stepping, the accepted
1047        // path should still integrate to the correct total noise.
1048        let expected_var = tf;
1049        for j in 0..n_interior {
1050            let empirical_var = sum_sq[j] / n_trajectories as f64;
1051            // Relatively loose tolerance due to finite sample size
1052            assert!(
1053                (empirical_var - expected_var).abs() < 0.06,
1054                "Component {} variance: expected ~{}, got {}",
1055                j,
1056                expected_var,
1057                empirical_var
1058            );
1059        }
1060    }
1061
1062    #[test]
1063    fn test_adaptive_options_builder() {
1064        let opts: SpdeOptions<f64> = SpdeOptions::default()
1065            .adaptive(true)
1066            .rtol(1e-5)
1067            .atol(1e-8)
1068            .dt_min(1e-12)
1069            .dt_max(0.05);
1070
1071        assert!(opts.adaptive);
1072        assert!((opts.rtol - 1e-5).abs() < 1e-15);
1073        assert!((opts.atol - 1e-8).abs() < 1e-18);
1074        assert!((opts.dt_min - 1e-12).abs() < 1e-22);
1075        assert!((opts.dt_max - 0.05).abs() < 1e-10);
1076    }
1077
1078    #[test]
1079    fn test_adaptive_finite_solution() {
1080        // Verify adaptive solver produces finite results even with larger noise
1081        let system = StochasticHeat {
1082            alpha: 0.05,
1083            sigma: 0.1,
1084        };
1085
1086        let grid = Grid1D::uniform(0.0, 1.0, 22); // 22 total -> 20 interior
1087        let u0: Vec<f64> = grid
1088            .interior_points()
1089            .iter()
1090            .map(|&x| (std::f64::consts::PI * x).sin())
1091            .collect();
1092
1093        let options = SpdeOptions::default()
1094            .adaptive(true)
1095            .dt(0.0005)
1096            .rtol(1e-3)
1097            .atol(1e-5)
1098            .dt_max(0.005)
1099            .seed(777);
1100
1101        let result = MolSdeSolver::solve(&system, 0.0, 0.1, &u0, &grid, &options).unwrap();
1102
1103        assert!(result.success);
1104        for y in &result.y {
1105            for &val in y {
1106                assert!(val.is_finite(), "Adaptive solution should be finite");
1107            }
1108        }
1109
1110        // Should have some rejected steps indicating adaptation
1111        // (Not required to pass but informative)
1112        assert!(result.stats.n_steps > 0, "Should have taken some steps");
1113    }
1114
1115    // ========================================================================
1116    // Colored / Trace-class noise integration tests
1117    // ========================================================================
1118
1119    struct ColoredHeat {
1120        alpha: f64,
1121        sigma: f64,
1122        correlation_length: f64,
1123    }
1124
1125    impl SpdeSystem<f64> for ColoredHeat {
1126        fn drift(&self, _t: f64, u: &[f64], du: &mut [f64], grid: &Grid1D<f64>) {
1127            let dx = grid.dx_uniform();
1128            let n = u.len();
1129            for i in 0..n {
1130                let u_left = if i == 0 { 0.0 } else { u[i - 1] };
1131                let u_right = if i == n - 1 { 0.0 } else { u[i + 1] };
1132                du[i] = self.alpha * (u_left - 2.0 * u[i] + u_right) / (dx * dx);
1133            }
1134        }
1135
1136        fn diffusion(&self, _t: f64, _u: &[f64], sigma: &mut [f64], _grid: &Grid1D<f64>) {
1137            for s in sigma.iter_mut() {
1138                *s = self.sigma;
1139            }
1140        }
1141
1142        fn noise_correlation(&self) -> crate::system::NoiseCorrelation<f64> {
1143            crate::system::NoiseCorrelation::Colored {
1144                correlation_length: self.correlation_length,
1145            }
1146        }
1147    }
1148
1149    #[test]
1150    fn test_spde_colored_noise_fixed_step() {
1151        // Stochastic heat with spatially correlated noise
1152        let system = ColoredHeat {
1153            alpha: 0.01,
1154            sigma: 0.1,
1155            correlation_length: 0.2,
1156        };
1157
1158        let grid = Grid1D::uniform(0.0, 1.0, 12); // 10 interior points
1159        let u0: Vec<f64> = grid
1160            .interior_points()
1161            .iter()
1162            .map(|&x| (std::f64::consts::PI * x).sin())
1163            .collect();
1164
1165        let options = SpdeOptions::default().dt(0.0001).n_output(20).seed(42);
1166
1167        let result = MolSdeSolver::solve(&system, 0.0, 0.01, &u0, &grid, &options).unwrap();
1168
1169        assert!(result.success);
1170        for y in &result.y {
1171            for &val in y {
1172                assert!(val.is_finite(), "Colored noise solution should be finite");
1173            }
1174        }
1175    }
1176
1177    #[test]
1178    fn test_spde_colored_noise_spatial_correlation() {
1179        // Verify that colored noise produces spatially correlated solutions.
1180        // Run an ensemble of pure-noise problems (zero drift) and check that
1181        // adjacent grid points are more correlated than distant ones.
1182        struct PureColoredNoise {
1183            correlation_length: f64,
1184        }
1185
1186        impl SpdeSystem<f64> for PureColoredNoise {
1187            fn drift(&self, _t: f64, _u: &[f64], du: &mut [f64], _grid: &Grid1D<f64>) {
1188                for d in du.iter_mut() {
1189                    *d = 0.0;
1190                }
1191            }
1192            fn diffusion(&self, _t: f64, _u: &[f64], sigma: &mut [f64], _grid: &Grid1D<f64>) {
1193                for s in sigma.iter_mut() {
1194                    *s = 1.0;
1195                }
1196            }
1197            fn noise_correlation(&self) -> crate::system::NoiseCorrelation<f64> {
1198                crate::system::NoiseCorrelation::Colored {
1199                    correlation_length: self.correlation_length,
1200                }
1201            }
1202        }
1203
1204        let system = PureColoredNoise {
1205            correlation_length: 0.3,
1206        };
1207        let grid = Grid1D::uniform(0.0, 1.0, 12); // 10 interior
1208        let n_interior = grid.len() - 2;
1209        let u0 = vec![0.0; n_interior];
1210
1211        let n_traj = 500;
1212        let mut finals = Vec::with_capacity(n_traj);
1213        for i in 0..n_traj {
1214            let options = SpdeOptions::default()
1215                .dt(0.001)
1216                .n_output(5)
1217                .seed(1000 + i as u64);
1218            let result = MolSdeSolver::solve(&system, 0.0, 0.01, &u0, &grid, &options).unwrap();
1219            finals.push(result.y_final().unwrap().to_vec());
1220        }
1221
1222        // Compute correlation between adjacent points (i=0, i=1)
1223        let mean_0: f64 = finals.iter().map(|f| f[0]).sum::<f64>() / n_traj as f64;
1224        let mean_1: f64 = finals.iter().map(|f| f[1]).sum::<f64>() / n_traj as f64;
1225        let var_0: f64 =
1226            finals.iter().map(|f| (f[0] - mean_0).powi(2)).sum::<f64>() / n_traj as f64;
1227        let var_1: f64 =
1228            finals.iter().map(|f| (f[1] - mean_1).powi(2)).sum::<f64>() / n_traj as f64;
1229        let cov_01: f64 = finals
1230            .iter()
1231            .map(|f| (f[0] - mean_0) * (f[1] - mean_1))
1232            .sum::<f64>()
1233            / n_traj as f64;
1234        let corr_nearby = cov_01 / (var_0.sqrt() * var_1.sqrt());
1235
1236        // Compute correlation between distant points (i=0, i=last)
1237        let last = n_interior - 1;
1238        let mean_far: f64 = finals.iter().map(|f| f[last]).sum::<f64>() / n_traj as f64;
1239        let var_far: f64 = finals
1240            .iter()
1241            .map(|f| (f[last] - mean_far).powi(2))
1242            .sum::<f64>()
1243            / n_traj as f64;
1244        let cov_far: f64 = finals
1245            .iter()
1246            .map(|f| (f[0] - mean_0) * (f[last] - mean_far))
1247            .sum::<f64>()
1248            / n_traj as f64;
1249        let corr_far = cov_far / (var_0.sqrt() * var_far.sqrt());
1250
1251        // Nearby correlation should be positive and larger than distant
1252        assert!(
1253            corr_nearby > 0.3,
1254            "Adjacent points should be positively correlated: {}",
1255            corr_nearby
1256        );
1257        assert!(
1258            corr_nearby > corr_far.abs(),
1259            "Adjacent correlation ({}) should exceed distant correlation ({})",
1260            corr_nearby,
1261            corr_far
1262        );
1263    }
1264
1265    #[test]
1266    fn test_spde_colored_noise_adaptive() {
1267        // Colored noise with adaptive time stepping
1268        let system = ColoredHeat {
1269            alpha: 0.01,
1270            sigma: 0.05,
1271            correlation_length: 0.2,
1272        };
1273
1274        let grid = Grid1D::uniform(0.0, 1.0, 12);
1275        let u0: Vec<f64> = grid
1276            .interior_points()
1277            .iter()
1278            .map(|&x| (std::f64::consts::PI * x).sin())
1279            .collect();
1280
1281        let options = SpdeOptions::default()
1282            .adaptive(true)
1283            .dt(0.001)
1284            .rtol(1e-3)
1285            .atol(1e-5)
1286            .dt_max(0.01)
1287            .seed(42);
1288
1289        let result = MolSdeSolver::solve(&system, 0.0, 0.05, &u0, &grid, &options).unwrap();
1290
1291        assert!(result.success);
1292        for y in &result.y {
1293            for &val in y {
1294                assert!(
1295                    val.is_finite(),
1296                    "Adaptive colored noise solution should be finite"
1297                );
1298            }
1299        }
1300    }
1301
1302    struct TraceClassHeat {
1303        alpha: f64,
1304        sigma: f64,
1305        n_modes: usize,
1306        decay_rate: f64,
1307    }
1308
1309    impl SpdeSystem<f64> for TraceClassHeat {
1310        fn drift(&self, _t: f64, u: &[f64], du: &mut [f64], grid: &Grid1D<f64>) {
1311            let dx = grid.dx_uniform();
1312            let n = u.len();
1313            for i in 0..n {
1314                let u_left = if i == 0 { 0.0 } else { u[i - 1] };
1315                let u_right = if i == n - 1 { 0.0 } else { u[i + 1] };
1316                du[i] = self.alpha * (u_left - 2.0 * u[i] + u_right) / (dx * dx);
1317            }
1318        }
1319
1320        fn diffusion(&self, _t: f64, _u: &[f64], sigma: &mut [f64], _grid: &Grid1D<f64>) {
1321            for s in sigma.iter_mut() {
1322                *s = self.sigma;
1323            }
1324        }
1325
1326        fn noise_correlation(&self) -> crate::system::NoiseCorrelation<f64> {
1327            crate::system::NoiseCorrelation::TraceClass {
1328                n_modes: self.n_modes,
1329                decay_rate: self.decay_rate,
1330            }
1331        }
1332    }
1333
1334    #[test]
1335    fn test_spde_trace_class_noise_fixed_step() {
1336        let system = TraceClassHeat {
1337            alpha: 0.01,
1338            sigma: 0.1,
1339            n_modes: 5,
1340            decay_rate: 2.0,
1341        };
1342
1343        let grid = Grid1D::uniform(0.0, 1.0, 22); // 20 interior
1344        let u0: Vec<f64> = grid
1345            .interior_points()
1346            .iter()
1347            .map(|&x| (std::f64::consts::PI * x).sin())
1348            .collect();
1349
1350        let options = SpdeOptions::default().dt(0.0001).n_output(20).seed(42);
1351
1352        let result = MolSdeSolver::solve(&system, 0.0, 0.01, &u0, &grid, &options).unwrap();
1353
1354        assert!(result.success);
1355        for y in &result.y {
1356            for &val in y {
1357                assert!(
1358                    val.is_finite(),
1359                    "Trace-class noise solution should be finite"
1360                );
1361            }
1362        }
1363    }
1364
1365    #[test]
1366    fn test_spde_trace_class_noise_smoother_than_white() {
1367        // Trace-class noise should produce smoother solutions than white noise
1368        // because higher modes are suppressed by eigenvalue decay.
1369        // Compare total variation of final solutions.
1370        struct PureNoise;
1371        impl SpdeSystem<f64> for PureNoise {
1372            fn drift(&self, _t: f64, _u: &[f64], du: &mut [f64], _grid: &Grid1D<f64>) {
1373                for d in du.iter_mut() {
1374                    *d = 0.0;
1375                }
1376            }
1377            fn diffusion(&self, _t: f64, _u: &[f64], sigma: &mut [f64], _grid: &Grid1D<f64>) {
1378                for s in sigma.iter_mut() {
1379                    *s = 1.0;
1380                }
1381            }
1382        }
1383
1384        struct PureTraceNoise;
1385        impl SpdeSystem<f64> for PureTraceNoise {
1386            fn drift(&self, _t: f64, _u: &[f64], du: &mut [f64], _grid: &Grid1D<f64>) {
1387                for d in du.iter_mut() {
1388                    *d = 0.0;
1389                }
1390            }
1391            fn diffusion(&self, _t: f64, _u: &[f64], sigma: &mut [f64], _grid: &Grid1D<f64>) {
1392                for s in sigma.iter_mut() {
1393                    *s = 1.0;
1394                }
1395            }
1396            fn noise_correlation(&self) -> crate::system::NoiseCorrelation<f64> {
1397                crate::system::NoiseCorrelation::TraceClass {
1398                    n_modes: 3,
1399                    decay_rate: 2.0,
1400                }
1401            }
1402        }
1403
1404        let grid = Grid1D::uniform(0.0, 1.0, 32); // 30 interior
1405        let n_interior = grid.len() - 2;
1406        let u0 = vec![0.0; n_interior];
1407
1408        let n_traj = 100;
1409        let mut white_tv_total = 0.0;
1410        let mut trace_tv_total = 0.0;
1411
1412        for i in 0..n_traj {
1413            let options = SpdeOptions::default()
1414                .dt(0.001)
1415                .n_output(5)
1416                .seed(2000 + i as u64);
1417
1418            let white_result =
1419                MolSdeSolver::solve(&PureNoise, 0.0, 0.01, &u0, &grid, &options).unwrap();
1420            let trace_result =
1421                MolSdeSolver::solve(&PureTraceNoise, 0.0, 0.01, &u0, &grid, &options).unwrap();
1422
1423            let white_y = white_result.y_final().unwrap();
1424            let trace_y = trace_result.y_final().unwrap();
1425
1426            // Total variation = sum |y[i+1] - y[i]|
1427            let white_tv: f64 = white_y.windows(2).map(|w| (w[1] - w[0]).abs()).sum();
1428            let trace_tv: f64 = trace_y.windows(2).map(|w| (w[1] - w[0]).abs()).sum();
1429
1430            white_tv_total += white_tv;
1431            trace_tv_total += trace_tv;
1432        }
1433
1434        // Trace-class noise should be smoother (lower total variation)
1435        assert!(
1436            trace_tv_total < white_tv_total,
1437            "Trace-class noise should be smoother: trace TV={}, white TV={}",
1438            trace_tv_total / n_traj as f64,
1439            white_tv_total / n_traj as f64
1440        );
1441    }
1442}