1use 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#[derive(Clone, Copy, Debug, Default)]
46pub enum SpdeMethod {
47 #[default]
49 EulerMaruyama,
50 Milstein,
52}
53
54#[derive(Clone, Debug)]
56pub struct SpdeOptions<S: Scalar> {
57 pub dt: S,
59 pub n_output: usize,
61 pub method: SpdeMethod,
63 pub seed: Option<u64>,
65 pub adaptive: bool,
67 pub rtol: S,
69 pub atol: S,
71 pub dt_min: S,
73 pub dt_max: S,
75}
76
77impl<S: Scalar> Default for SpdeOptions<S> {
78 fn default() -> Self {
79 Self {
80 dt: S::from_f64(0.001),
81 n_output: 100,
82 method: SpdeMethod::EulerMaruyama,
83 seed: None,
84 adaptive: false,
85 rtol: S::from_f64(1e-3),
86 atol: S::from_f64(1e-6),
87 dt_min: S::from_f64(1e-10),
88 dt_max: S::from_f64(0.1),
89 }
90 }
91}
92
93impl<S: Scalar> SpdeOptions<S> {
94 pub fn dt(mut self, dt: S) -> Self {
96 self.dt = dt;
97 self
98 }
99
100 pub fn n_output(mut self, n: usize) -> Self {
102 self.n_output = n;
103 self
104 }
105
106 pub fn method(mut self, method: SpdeMethod) -> Self {
108 self.method = method;
109 self
110 }
111
112 pub fn seed(mut self, seed: u64) -> Self {
114 self.seed = Some(seed);
115 self
116 }
117
118 pub fn adaptive(mut self, enable: bool) -> Self {
120 self.adaptive = enable;
121 self
122 }
123
124 pub fn rtol(mut self, tol: S) -> Self {
126 self.rtol = tol;
127 self
128 }
129
130 pub fn atol(mut self, tol: S) -> Self {
132 self.atol = tol;
133 self
134 }
135
136 pub fn dt_min(mut self, dt: S) -> Self {
138 self.dt_min = dt;
139 self
140 }
141
142 pub fn dt_max(mut self, dt: S) -> Self {
144 self.dt_max = dt;
145 self
146 }
147}
148
149#[derive(Clone, Debug, Default)]
151pub struct SpdeStats {
152 pub n_steps: usize,
154 pub n_drift: usize,
156 pub n_diffusion: usize,
158 pub n_reject: usize,
160}
161
162#[derive(Clone, Debug)]
164pub struct SpdeResult<S: Scalar> {
165 pub t: Vec<S>,
167 pub y: Vec<Vec<S>>,
169 pub grid: Grid1D<S>,
171 pub stats: SpdeStats,
173 pub success: bool,
175}
176
177impl<S: Scalar> SpdeResult<S> {
178 pub fn y_final(&self) -> Option<&[S]> {
180 self.y.last().map(|v| v.as_slice())
181 }
182
183 pub fn t_final(&self) -> Option<S> {
185 self.t.last().copied()
186 }
187
188 pub fn y_at(&self, i: usize) -> &[S] {
190 &self.y[i]
191 }
192}
193
194pub trait SpdeSolver<S: Scalar> {
196 fn solve<Sys: SpdeSystem<S> + Sync>(
198 system: &Sys,
199 t0: S,
200 tf: S,
201 u0: &[S],
202 grid: &Grid1D<S>,
203 options: &SpdeOptions<S>,
204 ) -> Result<SpdeResult<S>, String>;
205}
206
207enum MolNoiseConfig<S: Scalar> {
209 White,
211 Colored { cholesky: Vec<S>, dim: usize },
214 TraceClass {
217 sqrt_eigenvalues: Vec<S>,
219 basis: Vec<S>,
221 n_modes: usize,
222 dim: usize,
223 },
224}
225
226impl<S: Scalar> MolNoiseConfig<S> {
227 fn from_system<Sys: SpdeSystem<S>>(system: &Sys, grid: &Grid1D<S>) -> Result<Self, String> {
229 match system.noise_correlation() {
230 NoiseCorrelation::White => Ok(MolNoiseConfig::White),
231 NoiseCorrelation::Colored { correlation_length } => {
232 let dim = grid.len() - 2; let dx = grid.dx_uniform();
234
235 let mut c = vec![S::ZERO; dim * dim];
237 for i in 0..dim {
238 for j in 0..dim {
239 let dist = dx * S::from_usize(if i > j { i - j } else { j - i });
240 c[i * dim + j] = (-dist / correlation_length).exp();
241 }
242 }
243
244 let cholesky = cholesky_decompose(&c, dim)?;
245 Ok(MolNoiseConfig::Colored { cholesky, dim })
246 }
247 NoiseCorrelation::TraceClass {
248 n_modes,
249 decay_rate,
250 } => {
251 let dim = grid.len() - 2;
252 let interior = grid.interior_points();
253 let (x_min, x_max) = grid.bounds();
254 let domain_length = x_max - x_min;
255
256 let sqrt_eigenvalues: Vec<S> = (1..=n_modes)
258 .map(|k| S::from_usize(k).powf(-decay_rate).sqrt())
259 .collect();
260
261 let norm = (S::TWO / domain_length).sqrt();
263 let pi = S::PI;
264 let mut basis = vec![S::ZERO; dim * n_modes];
265 for i in 0..dim {
266 let x = interior[i] - x_min;
267 for k in 0..n_modes {
268 let mode = S::from_usize(k + 1);
269 basis[i * n_modes + k] = norm * (mode * pi * x / domain_length).sin();
270 }
271 }
272
273 Ok(MolNoiseConfig::TraceClass {
274 sqrt_eigenvalues,
275 basis,
276 n_modes,
277 dim,
278 })
279 }
280 }
281 }
282
283 fn n_wiener(&self, dim: usize) -> usize {
285 match self {
286 MolNoiseConfig::White => dim,
287 MolNoiseConfig::Colored { dim: d, .. } => *d,
288 MolNoiseConfig::TraceClass { n_modes, .. } => *n_modes,
289 }
290 }
291}
292
293struct MolSdeWrapper<'a, S: Scalar, Sys: SpdeSystem<S> + Sync> {
295 spde: &'a Sys,
296 grid: &'a Grid1D<S>,
297 noise_config: MolNoiseConfig<S>,
298}
299
300impl<'a, S: Scalar, Sys: SpdeSystem<S> + Sync> SdeSystem<S> for MolSdeWrapper<'a, S, Sys> {
301 fn dim(&self) -> usize {
302 self.grid.len() - 2 }
304
305 fn noise_type(&self) -> NoiseType {
306 match &self.noise_config {
307 MolNoiseConfig::White => NoiseType::Diagonal,
308 MolNoiseConfig::Colored { dim, .. } => NoiseType::General { n_wiener: *dim },
309 MolNoiseConfig::TraceClass { n_modes, .. } => NoiseType::General { n_wiener: *n_modes },
310 }
311 }
312
313 fn drift(&self, t: S, y: &[S], f: &mut [S]) {
314 self.spde.drift(t, y, f, self.grid);
315 }
316
317 fn diffusion(&self, t: S, y: &[S], g: &mut [S]) {
318 match &self.noise_config {
319 MolNoiseConfig::White => {
320 self.spde.diffusion(t, y, g, self.grid);
322 }
323 MolNoiseConfig::Colored { cholesky, dim } => {
324 let mut sigma = vec![S::ZERO; *dim];
327 self.spde.diffusion(t, y, &mut sigma, self.grid);
328
329 for i in 0..*dim {
331 for j in 0..=i {
332 g[i * dim + j] = sigma[i] * cholesky[i * dim + j];
333 }
334 for j in (i + 1)..*dim {
336 g[i * dim + j] = S::ZERO;
337 }
338 }
339 }
340 MolNoiseConfig::TraceClass {
341 sqrt_eigenvalues,
342 basis,
343 n_modes,
344 dim,
345 } => {
346 let mut sigma = vec![S::ZERO; *dim];
348 self.spde.diffusion(t, y, &mut sigma, self.grid);
349
350 for i in 0..*dim {
351 for k in 0..*n_modes {
352 g[i * n_modes + k] =
353 sigma[i] * sqrt_eigenvalues[k] * basis[i * n_modes + k];
354 }
355 }
356 }
357 }
358 }
359}
360
361unsafe impl<'a, S: Scalar, Sys: SpdeSystem<S> + Sync> Sync for MolSdeWrapper<'a, S, Sys> {}
369
370pub struct MolSdeSolver;
375
376impl<S: Scalar> SpdeSolver<S> for MolSdeSolver {
377 fn solve<Sys: SpdeSystem<S> + Sync>(
378 system: &Sys,
379 t0: S,
380 tf: S,
381 u0: &[S],
382 grid: &Grid1D<S>,
383 options: &SpdeOptions<S>,
384 ) -> Result<SpdeResult<S>, String> {
385 let n_interior = grid.len() - 2;
386 if u0.len() != n_interior {
387 return Err(format!(
388 "Initial condition length {} doesn't match interior points {}",
389 u0.len(),
390 n_interior
391 ));
392 }
393
394 if options.adaptive {
396 Self::solve_adaptive(system, t0, tf, u0, grid, options)
397 } else {
398 Self::solve_fixed(system, t0, tf, u0, grid, options)
399 }
400 }
401}
402
403impl MolSdeSolver {
404 fn solve_fixed<S: Scalar, Sys: SpdeSystem<S> + Sync>(
406 system: &Sys,
407 t0: S,
408 tf: S,
409 u0: &[S],
410 grid: &Grid1D<S>,
411 options: &SpdeOptions<S>,
412 ) -> Result<SpdeResult<S>, String> {
413 let noise_config = MolNoiseConfig::from_system(system, grid)?;
415 let sde_system = MolSdeWrapper {
416 spde: system,
417 grid,
418 noise_config,
419 };
420
421 let sde_options = SdeOptions::default().dt(options.dt);
423
424 let sde_result = match options.method {
426 SpdeMethod::EulerMaruyama => {
427 EulerMaruyama::solve(&sde_system, t0, tf, u0, &sde_options, options.seed)
428 }
429 SpdeMethod::Milstein => {
430 if system.is_additive() {
431 EulerMaruyama::solve(&sde_system, t0, tf, u0, &sde_options, options.seed)
433 } else {
434 Milstein::solve(&sde_system, t0, tf, u0, &sde_options, options.seed)
435 }
436 }
437 }?;
438
439 let n_steps = sde_result.t.len();
441 let stride = (n_steps / options.n_output).max(1);
442
443 let mut t_out = Vec::new();
444 let mut y_out = Vec::new();
445
446 for i in (0..n_steps).step_by(stride) {
447 t_out.push(sde_result.t[i]);
448 y_out.push(sde_result.y_at(i).to_vec());
449 }
450
451 if t_out.last() != sde_result.t.last() {
453 t_out.push(*sde_result.t.last().unwrap());
454 y_out.push(sde_result.y_final().unwrap().to_vec());
455 }
456
457 Ok(SpdeResult {
458 t: t_out,
459 y: y_out,
460 grid: grid.clone(),
461 stats: SpdeStats {
462 n_steps: sde_result.stats.n_accept + sde_result.stats.n_reject,
463 n_drift: sde_result.stats.n_drift,
464 n_diffusion: sde_result.stats.n_diffusion,
465 n_reject: sde_result.stats.n_reject,
466 },
467 success: sde_result.success,
468 })
469 }
470
471 fn solve_adaptive<S: Scalar, Sys: SpdeSystem<S> + Sync>(
473 system: &Sys,
474 t0: S,
475 tf: S,
476 u0: &[S],
477 grid: &Grid1D<S>,
478 options: &SpdeOptions<S>,
479 ) -> Result<SpdeResult<S>, String> {
480 let dim = u0.len();
481 let mut t = t0;
482 let mut u = u0.to_vec();
483 let mut dt = options.dt;
484
485 let noise_config = MolNoiseConfig::from_system(system, grid)?;
487 let n_wiener = noise_config.n_wiener(dim);
488
489 let mut wiener = create_wiener(n_wiener, options.seed);
491
492 let mut t_out = vec![t0];
494 let mut y_out = vec![u0.to_vec()];
495 let mut stats = SpdeStats::default();
496
497 let mut drift = vec![S::ZERO; dim];
499 let mut diffusion = vec![S::ZERO; dim];
500 let mut u1 = vec![S::ZERO; dim]; let mut u_half = vec![S::ZERO; dim]; let mut u2 = vec![S::ZERO; dim]; let max_steps = 1_000_000;
505 let mut step_count = 0;
506
507 while t < tf && step_count < max_steps {
508 let h = dt.min(tf - t).min(options.dt_max).max(options.dt_min);
510
511 let half_h = h / S::from_f64(2.0);
514 let raw_dw1 = wiener.increment(half_h);
515 let raw_dw2 = wiener.increment(half_h);
516
517 let noise1 = apply_noise_correlation(&noise_config, &raw_dw1.dw, dim);
521 let noise2 = apply_noise_correlation(&noise_config, &raw_dw2.dw, dim);
522
523 let mut drift_orig = vec![S::ZERO; dim];
525 let mut diffusion_orig = vec![S::ZERO; dim];
526 system.drift(t, &u, &mut drift_orig, grid);
527 system.diffusion(t, &u, &mut diffusion_orig, grid);
528 stats.n_drift += 1;
529 stats.n_diffusion += 1;
530
531 for i in 0..dim {
533 u1[i] = u[i] + drift_orig[i] * h + diffusion_orig[i] * (noise1[i] + noise2[i]);
534 }
535
536 for i in 0..dim {
538 u_half[i] = u[i] + drift_orig[i] * half_h + diffusion_orig[i] * noise1[i];
539 }
540
541 system.drift(t + half_h, &u_half, &mut drift, grid);
543 system.diffusion(t + half_h, &u_half, &mut diffusion, grid);
544 stats.n_drift += 1;
545 stats.n_diffusion += 1;
546
547 for i in 0..dim {
549 u2[i] = u_half[i] + drift[i] * half_h + diffusion[i] * noise2[i];
550 }
551
552 let mut err_max = S::ZERO;
554 let mut scale_max = S::ZERO;
555 for i in 0..dim {
556 let err_i = (u1[i] - u2[i]).abs();
557 let scale_i = options.atol + options.rtol * u[i].abs().max(u2[i].abs());
558 let ratio = err_i / scale_i;
559 if ratio > err_max / scale_max.max(S::from_f64(1e-15)) {
560 err_max = err_i;
561 scale_max = scale_i;
562 }
563 }
564
565 let err_ratio = err_max / scale_max.max(S::from_f64(1e-15));
566
567 if err_ratio <= S::ONE || h <= options.dt_min {
568 t = t + h;
570 u.copy_from_slice(&u2);
571 stats.n_steps += 1;
572
573 t_out.push(t);
575 y_out.push(u.clone());
576
577 if err_ratio < S::from_f64(0.5) {
580 let factor = (S::ONE / err_ratio.max(S::from_f64(1e-10)))
581 .powf(S::from_f64(1.0 / 3.0))
582 .min(S::from_f64(2.0));
583 dt = (h * factor).min(options.dt_max);
584 }
585 } else {
586 stats.n_reject += 1;
588 let factor = (S::ONE / err_ratio.max(S::from_f64(1e-10)))
589 .powf(S::from_f64(1.0 / 3.0))
590 .max(S::from_f64(0.25))
591 .min(S::from_f64(0.9));
592 dt = (h * factor).max(options.dt_min);
593 }
594
595 step_count += 1;
596 }
597
598 if step_count >= max_steps && t < tf {
599 return Err(format!(
600 "Maximum steps ({}) exceeded at t = {}",
601 max_steps,
602 t.to_f64()
603 ));
604 }
605
606 let n_total = t_out.len();
608 if n_total > options.n_output {
609 let stride = (n_total / options.n_output).max(1);
610 let mut t_sub = Vec::new();
611 let mut y_sub = Vec::new();
612
613 for i in (0..n_total).step_by(stride) {
614 t_sub.push(t_out[i]);
615 y_sub.push(y_out[i].clone());
616 }
617
618 if t_sub.last() != t_out.last() {
620 t_sub.push(*t_out.last().unwrap());
621 y_sub.push(y_out.last().unwrap().clone());
622 }
623
624 t_out = t_sub;
625 y_out = y_sub;
626 }
627
628 Ok(SpdeResult {
629 t: t_out,
630 y: y_out,
631 grid: grid.clone(),
632 stats,
633 success: true,
634 })
635 }
636}
637
638fn apply_noise_correlation<S: Scalar>(
644 config: &MolNoiseConfig<S>,
645 raw_dw: &[S],
646 dim: usize,
647) -> Vec<S> {
648 match config {
649 MolNoiseConfig::White => {
650 raw_dw[..dim].to_vec()
652 }
653 MolNoiseConfig::Colored { cholesky, dim: n } => {
654 let mut result = vec![S::ZERO; *n];
656 for i in 0..*n {
657 for j in 0..=i {
658 result[i] = result[i] + cholesky[i * n + j] * raw_dw[j];
659 }
660 }
661 result
662 }
663 MolNoiseConfig::TraceClass {
664 sqrt_eigenvalues,
665 basis,
666 n_modes,
667 dim: n,
668 } => {
669 let mut result = vec![S::ZERO; *n];
671 for i in 0..*n {
672 for k in 0..*n_modes {
673 result[i] =
674 result[i] + sqrt_eigenvalues[k] * basis[i * n_modes + k] * raw_dw[k];
675 }
676 }
677 result
678 }
679 }
680}
681
682pub struct SpdeEnsemble;
684
685impl SpdeEnsemble {
686 pub fn solve<S: Scalar, Sys: SpdeSystem<S> + Sync>(
688 system: &Sys,
689 t0: S,
690 tf: S,
691 u0: &[S],
692 grid: &Grid1D<S>,
693 options: &SpdeOptions<S>,
694 n_trajectories: usize,
695 ) -> Result<Vec<SpdeResult<S>>, String> {
696 let mut results = Vec::with_capacity(n_trajectories);
697 let base_seed = options.seed.unwrap_or(42);
698
699 for i in 0..n_trajectories {
700 let mut opts = options.clone();
701 opts.seed = Some(base_seed + i as u64);
702 let result = MolSdeSolver::solve(system, t0, tf, u0, grid, &opts)?;
703 results.push(result);
704 }
705
706 Ok(results)
707 }
708
709 pub fn mean<S: Scalar>(results: &[SpdeResult<S>]) -> Vec<Vec<S>> {
711 if results.is_empty() {
712 return Vec::new();
713 }
714
715 let n_times = results[0].t.len();
716 let n_points = results[0].y[0].len();
717 let n_traj = results.len();
718
719 let mut mean = vec![vec![S::ZERO; n_points]; n_times];
720
721 for result in results {
722 for (i, y) in result.y.iter().enumerate() {
723 for (j, &val) in y.iter().enumerate() {
724 mean[i][j] = mean[i][j] + val;
725 }
726 }
727 }
728
729 let n_traj_s = S::from_usize(n_traj);
730 for m in mean.iter_mut() {
731 for val in m.iter_mut() {
732 *val = *val / n_traj_s;
733 }
734 }
735
736 mean
737 }
738
739 pub fn std<S: Scalar>(results: &[SpdeResult<S>], mean: &[Vec<S>]) -> Vec<Vec<S>> {
741 if results.is_empty() {
742 return Vec::new();
743 }
744
745 let n_times = results[0].t.len();
746 let n_points = results[0].y[0].len();
747 let n_traj = results.len();
748
749 let mut variance = vec![vec![S::ZERO; n_points]; n_times];
750
751 for result in results {
752 for (i, y) in result.y.iter().enumerate() {
753 for (j, &val) in y.iter().enumerate() {
754 let diff = val - mean[i][j];
755 variance[i][j] = variance[i][j] + diff * diff;
756 }
757 }
758 }
759
760 let n_traj_s = S::from_usize(n_traj - 1); for v in variance.iter_mut() {
762 for val in v.iter_mut() {
763 *val = (*val / n_traj_s).sqrt();
764 }
765 }
766
767 variance
768 }
769}
770
771#[cfg(test)]
772mod tests {
773 use super::*;
774
775 struct StochasticHeat {
776 alpha: f64,
777 sigma: f64,
778 }
779
780 impl SpdeSystem<f64> for StochasticHeat {
781 fn drift(&self, _t: f64, u: &[f64], du: &mut [f64], grid: &Grid1D<f64>) {
782 let dx = grid.dx_uniform();
783 let n = u.len();
784
785 for i in 0..n {
786 let u_left = if i == 0 { 0.0 } else { u[i - 1] };
787 let u_right = if i == n - 1 { 0.0 } else { u[i + 1] };
788 du[i] = self.alpha * (u_left - 2.0 * u[i] + u_right) / (dx * dx);
789 }
790 }
791
792 fn diffusion(&self, _t: f64, _u: &[f64], sigma: &mut [f64], _grid: &Grid1D<f64>) {
793 for s in sigma.iter_mut() {
794 *s = self.sigma;
795 }
796 }
797 }
798
799 #[test]
800 fn test_stochastic_heat_basic() {
801 let system = StochasticHeat {
802 alpha: 0.01,
803 sigma: 0.1,
804 };
805
806 let grid = Grid1D::uniform(0.0, 1.0, 21);
807 let u0: Vec<f64> = grid
808 .interior_points()
809 .iter()
810 .map(|&x| (std::f64::consts::PI * x).sin())
811 .collect();
812
813 let options = SpdeOptions::default().dt(0.0001).n_output(50).seed(12345);
814
815 let result = MolSdeSolver::solve(&system, 0.0, 0.1, &u0, &grid, &options).unwrap();
816
817 assert!(result.success);
818 assert!(!result.t.is_empty());
819 assert!(!result.y.is_empty());
820
821 for y in &result.y {
823 for &val in y {
824 assert!(val.is_finite(), "Solution should be finite");
825 }
826 }
827 }
828
829 #[test]
830 fn test_spde_reproducibility() {
831 let system = StochasticHeat {
832 alpha: 0.01,
833 sigma: 0.1,
834 };
835
836 let grid = Grid1D::uniform(0.0, 1.0, 11);
837 let u0: Vec<f64> = grid
838 .interior_points()
839 .iter()
840 .map(|&x| (std::f64::consts::PI * x).sin())
841 .collect();
842
843 let options = SpdeOptions::default().dt(0.001).n_output(10).seed(42);
844
845 let result1 = MolSdeSolver::solve(&system, 0.0, 0.01, &u0, &grid, &options).unwrap();
846 let result2 = MolSdeSolver::solve(&system, 0.0, 0.01, &u0, &grid, &options).unwrap();
847
848 let y1 = result1.y_final().unwrap();
849 let y2 = result2.y_final().unwrap();
850
851 for (a, b) in y1.iter().zip(y2.iter()) {
852 assert!(
853 (a - b).abs() < 1e-10,
854 "Results should be reproducible with same seed"
855 );
856 }
857 }
858
859 #[test]
860 fn test_ensemble_statistics() {
861 let system = StochasticHeat {
862 alpha: 0.01,
863 sigma: 0.05,
864 };
865
866 let grid = Grid1D::uniform(0.0, 1.0, 11);
867 let u0: Vec<f64> = grid
868 .interior_points()
869 .iter()
870 .map(|&x| (std::f64::consts::PI * x).sin())
871 .collect();
872
873 let options = SpdeOptions::default().dt(0.001).n_output(10).seed(100);
874
875 let results = SpdeEnsemble::solve(&system, 0.0, 0.01, &u0, &grid, &options, 10).unwrap();
876
877 assert_eq!(results.len(), 10);
878
879 let mean = SpdeEnsemble::mean(&results);
880 let std = SpdeEnsemble::std(&results, &mean);
881
882 for s in &std[std.len() - 1] {
885 assert!(*s >= 0.0, "Std should be non-negative");
886 }
887 }
888
889 #[test]
890 fn test_zero_noise() {
891 let system = StochasticHeat {
893 alpha: 0.1,
894 sigma: 0.0,
895 };
896
897 let grid = Grid1D::uniform(0.0, 1.0, 11);
898 let u0: Vec<f64> = grid
899 .interior_points()
900 .iter()
901 .map(|&x| (std::f64::consts::PI * x).sin())
902 .collect();
903
904 let options = SpdeOptions::default().dt(0.0001).n_output(10).seed(42);
905
906 let result1 = MolSdeSolver::solve(&system, 0.0, 0.01, &u0, &grid, &options).unwrap();
907
908 let options2 = options.seed(999);
910 let result2 = MolSdeSolver::solve(&system, 0.0, 0.01, &u0, &grid, &options2).unwrap();
911
912 let y1 = result1.y_final().unwrap();
913 let y2 = result2.y_final().unwrap();
914
915 for (a, b) in y1.iter().zip(y2.iter()) {
916 assert!(
917 (a - b).abs() < 1e-8,
918 "Zero noise should give deterministic results"
919 );
920 }
921 }
922
923 #[test]
924 fn test_options_builder() {
925 let opts: SpdeOptions<f64> = SpdeOptions::default()
926 .dt(0.01)
927 .n_output(200)
928 .method(SpdeMethod::Milstein)
929 .seed(123);
930
931 assert!((opts.dt - 0.01).abs() < 1e-10);
932 assert_eq!(opts.n_output, 200);
933 assert!(matches!(opts.method, SpdeMethod::Milstein));
934 assert_eq!(opts.seed, Some(123));
935 }
936
937 #[test]
938 fn test_spde_adaptive_dt() {
939 let spde = StochasticHeat {
941 alpha: 0.1,
942 sigma: 0.01,
943 };
944
945 let grid = Grid1D::uniform(0.0, 1.0, 52); let u0: Vec<f64> = grid
947 .interior_points()
948 .iter()
949 .map(|&x| (std::f64::consts::PI * x).sin())
950 .collect();
951
952 let options = SpdeOptions::default()
953 .adaptive(true)
954 .dt(0.001)
955 .rtol(1e-4)
956 .atol(1e-6)
957 .dt_max(0.01)
958 .seed(42);
959
960 let result = MolSdeSolver::solve(&spde, 0.0, 0.5, &u0, &grid, &options).unwrap();
961
962 assert!(result.success, "Adaptive solve should succeed");
963 assert!(!result.t.is_empty(), "Should have output time points");
964
965 for y in &result.y {
967 for &val in y {
968 assert!(val.is_finite(), "Solution should be finite, got {}", val);
969 }
970 }
971
972 let dts: Vec<_> = result.t.windows(2).map(|w| w[1] - w[0]).collect();
974 let first_dt = dts.first().copied().unwrap_or(0.0);
975 let has_variation = dts.iter().any(|&dt| (dt - first_dt).abs() > 1e-10);
976
977 assert!(
978 has_variation,
979 "Time steps should vary with adaptive stepping"
980 );
981 }
982
983 #[test]
984 fn test_adaptive_noise_consistency() {
985 struct PureNoise;
992
993 impl SpdeSystem<f64> for PureNoise {
994 fn drift(&self, _t: f64, _u: &[f64], du: &mut [f64], _grid: &Grid1D<f64>) {
995 for d in du.iter_mut() {
996 *d = 0.0;
997 }
998 }
999
1000 fn diffusion(&self, _t: f64, _u: &[f64], sigma: &mut [f64], _grid: &Grid1D<f64>) {
1001 for s in sigma.iter_mut() {
1002 *s = 1.0;
1003 }
1004 }
1005 }
1006
1007 let n_trajectories = 500;
1008 let grid = Grid1D::uniform(0.0, 1.0, 12); let n_interior = grid.len() - 2;
1010 let u0 = vec![0.0; n_interior];
1011 let tf = 0.1;
1012
1013 let mut sum_sq = vec![0.0; n_interior];
1014 for i in 0..n_trajectories {
1015 let options = SpdeOptions::default()
1016 .adaptive(true)
1017 .dt(0.01)
1018 .rtol(1e-2)
1019 .atol(1e-4)
1020 .dt_min(1e-6)
1021 .dt_max(0.05)
1022 .seed(1000 + i as u64);
1023
1024 let result = MolSdeSolver::solve(&PureNoise, 0.0, tf, &u0, &grid, &options).unwrap();
1025 let y_final = result.y_final().unwrap();
1026 for j in 0..n_interior {
1027 sum_sq[j] += y_final[j] * y_final[j];
1028 }
1029 }
1030
1031 let expected_var = tf;
1035 for j in 0..n_interior {
1036 let empirical_var = sum_sq[j] / n_trajectories as f64;
1037 assert!(
1039 (empirical_var - expected_var).abs() < 0.06,
1040 "Component {} variance: expected ~{}, got {}",
1041 j,
1042 expected_var,
1043 empirical_var
1044 );
1045 }
1046 }
1047
1048 #[test]
1049 fn test_adaptive_options_builder() {
1050 let opts: SpdeOptions<f64> = SpdeOptions::default()
1051 .adaptive(true)
1052 .rtol(1e-5)
1053 .atol(1e-8)
1054 .dt_min(1e-12)
1055 .dt_max(0.05);
1056
1057 assert!(opts.adaptive);
1058 assert!((opts.rtol - 1e-5).abs() < 1e-15);
1059 assert!((opts.atol - 1e-8).abs() < 1e-18);
1060 assert!((opts.dt_min - 1e-12).abs() < 1e-22);
1061 assert!((opts.dt_max - 0.05).abs() < 1e-10);
1062 }
1063
1064 #[test]
1065 fn test_adaptive_finite_solution() {
1066 let system = StochasticHeat {
1068 alpha: 0.05,
1069 sigma: 0.1,
1070 };
1071
1072 let grid = Grid1D::uniform(0.0, 1.0, 22); let u0: Vec<f64> = grid
1074 .interior_points()
1075 .iter()
1076 .map(|&x| (std::f64::consts::PI * x).sin())
1077 .collect();
1078
1079 let options = SpdeOptions::default()
1080 .adaptive(true)
1081 .dt(0.0005)
1082 .rtol(1e-3)
1083 .atol(1e-5)
1084 .dt_max(0.005)
1085 .seed(777);
1086
1087 let result = MolSdeSolver::solve(&system, 0.0, 0.1, &u0, &grid, &options).unwrap();
1088
1089 assert!(result.success);
1090 for y in &result.y {
1091 for &val in y {
1092 assert!(val.is_finite(), "Adaptive solution should be finite");
1093 }
1094 }
1095
1096 assert!(result.stats.n_steps > 0, "Should have taken some steps");
1099 }
1100
1101 struct ColoredHeat {
1106 alpha: f64,
1107 sigma: f64,
1108 correlation_length: f64,
1109 }
1110
1111 impl SpdeSystem<f64> for ColoredHeat {
1112 fn drift(&self, _t: f64, u: &[f64], du: &mut [f64], grid: &Grid1D<f64>) {
1113 let dx = grid.dx_uniform();
1114 let n = u.len();
1115 for i in 0..n {
1116 let u_left = if i == 0 { 0.0 } else { u[i - 1] };
1117 let u_right = if i == n - 1 { 0.0 } else { u[i + 1] };
1118 du[i] = self.alpha * (u_left - 2.0 * u[i] + u_right) / (dx * dx);
1119 }
1120 }
1121
1122 fn diffusion(&self, _t: f64, _u: &[f64], sigma: &mut [f64], _grid: &Grid1D<f64>) {
1123 for s in sigma.iter_mut() {
1124 *s = self.sigma;
1125 }
1126 }
1127
1128 fn noise_correlation(&self) -> crate::system::NoiseCorrelation<f64> {
1129 crate::system::NoiseCorrelation::Colored {
1130 correlation_length: self.correlation_length,
1131 }
1132 }
1133 }
1134
1135 #[test]
1136 fn test_spde_colored_noise_fixed_step() {
1137 let system = ColoredHeat {
1139 alpha: 0.01,
1140 sigma: 0.1,
1141 correlation_length: 0.2,
1142 };
1143
1144 let grid = Grid1D::uniform(0.0, 1.0, 12); let u0: Vec<f64> = grid
1146 .interior_points()
1147 .iter()
1148 .map(|&x| (std::f64::consts::PI * x).sin())
1149 .collect();
1150
1151 let options = SpdeOptions::default().dt(0.0001).n_output(20).seed(42);
1152
1153 let result = MolSdeSolver::solve(&system, 0.0, 0.01, &u0, &grid, &options).unwrap();
1154
1155 assert!(result.success);
1156 for y in &result.y {
1157 for &val in y {
1158 assert!(val.is_finite(), "Colored noise solution should be finite");
1159 }
1160 }
1161 }
1162
1163 #[test]
1164 fn test_spde_colored_noise_spatial_correlation() {
1165 struct PureColoredNoise {
1169 correlation_length: f64,
1170 }
1171
1172 impl SpdeSystem<f64> for PureColoredNoise {
1173 fn drift(&self, _t: f64, _u: &[f64], du: &mut [f64], _grid: &Grid1D<f64>) {
1174 for d in du.iter_mut() {
1175 *d = 0.0;
1176 }
1177 }
1178 fn diffusion(&self, _t: f64, _u: &[f64], sigma: &mut [f64], _grid: &Grid1D<f64>) {
1179 for s in sigma.iter_mut() {
1180 *s = 1.0;
1181 }
1182 }
1183 fn noise_correlation(&self) -> crate::system::NoiseCorrelation<f64> {
1184 crate::system::NoiseCorrelation::Colored {
1185 correlation_length: self.correlation_length,
1186 }
1187 }
1188 }
1189
1190 let system = PureColoredNoise {
1191 correlation_length: 0.3,
1192 };
1193 let grid = Grid1D::uniform(0.0, 1.0, 12); let n_interior = grid.len() - 2;
1195 let u0 = vec![0.0; n_interior];
1196
1197 let n_traj = 500;
1198 let mut finals = Vec::with_capacity(n_traj);
1199 for i in 0..n_traj {
1200 let options = SpdeOptions::default()
1201 .dt(0.001)
1202 .n_output(5)
1203 .seed(1000 + i as u64);
1204 let result = MolSdeSolver::solve(&system, 0.0, 0.01, &u0, &grid, &options).unwrap();
1205 finals.push(result.y_final().unwrap().to_vec());
1206 }
1207
1208 let mean_0: f64 = finals.iter().map(|f| f[0]).sum::<f64>() / n_traj as f64;
1210 let mean_1: f64 = finals.iter().map(|f| f[1]).sum::<f64>() / n_traj as f64;
1211 let var_0: f64 =
1212 finals.iter().map(|f| (f[0] - mean_0).powi(2)).sum::<f64>() / n_traj as f64;
1213 let var_1: f64 =
1214 finals.iter().map(|f| (f[1] - mean_1).powi(2)).sum::<f64>() / n_traj as f64;
1215 let cov_01: f64 = finals
1216 .iter()
1217 .map(|f| (f[0] - mean_0) * (f[1] - mean_1))
1218 .sum::<f64>()
1219 / n_traj as f64;
1220 let corr_nearby = cov_01 / (var_0.sqrt() * var_1.sqrt());
1221
1222 let last = n_interior - 1;
1224 let mean_far: f64 = finals.iter().map(|f| f[last]).sum::<f64>() / n_traj as f64;
1225 let var_far: f64 = finals
1226 .iter()
1227 .map(|f| (f[last] - mean_far).powi(2))
1228 .sum::<f64>()
1229 / n_traj as f64;
1230 let cov_far: f64 = finals
1231 .iter()
1232 .map(|f| (f[0] - mean_0) * (f[last] - mean_far))
1233 .sum::<f64>()
1234 / n_traj as f64;
1235 let corr_far = cov_far / (var_0.sqrt() * var_far.sqrt());
1236
1237 assert!(
1239 corr_nearby > 0.3,
1240 "Adjacent points should be positively correlated: {}",
1241 corr_nearby
1242 );
1243 assert!(
1244 corr_nearby > corr_far.abs(),
1245 "Adjacent correlation ({}) should exceed distant correlation ({})",
1246 corr_nearby,
1247 corr_far
1248 );
1249 }
1250
1251 #[test]
1252 fn test_spde_colored_noise_adaptive() {
1253 let system = ColoredHeat {
1255 alpha: 0.01,
1256 sigma: 0.05,
1257 correlation_length: 0.2,
1258 };
1259
1260 let grid = Grid1D::uniform(0.0, 1.0, 12);
1261 let u0: Vec<f64> = grid
1262 .interior_points()
1263 .iter()
1264 .map(|&x| (std::f64::consts::PI * x).sin())
1265 .collect();
1266
1267 let options = SpdeOptions::default()
1268 .adaptive(true)
1269 .dt(0.001)
1270 .rtol(1e-3)
1271 .atol(1e-5)
1272 .dt_max(0.01)
1273 .seed(42);
1274
1275 let result = MolSdeSolver::solve(&system, 0.0, 0.05, &u0, &grid, &options).unwrap();
1276
1277 assert!(result.success);
1278 for y in &result.y {
1279 for &val in y {
1280 assert!(
1281 val.is_finite(),
1282 "Adaptive colored noise solution should be finite"
1283 );
1284 }
1285 }
1286 }
1287
1288 struct TraceClassHeat {
1289 alpha: f64,
1290 sigma: f64,
1291 n_modes: usize,
1292 decay_rate: f64,
1293 }
1294
1295 impl SpdeSystem<f64> for TraceClassHeat {
1296 fn drift(&self, _t: f64, u: &[f64], du: &mut [f64], grid: &Grid1D<f64>) {
1297 let dx = grid.dx_uniform();
1298 let n = u.len();
1299 for i in 0..n {
1300 let u_left = if i == 0 { 0.0 } else { u[i - 1] };
1301 let u_right = if i == n - 1 { 0.0 } else { u[i + 1] };
1302 du[i] = self.alpha * (u_left - 2.0 * u[i] + u_right) / (dx * dx);
1303 }
1304 }
1305
1306 fn diffusion(&self, _t: f64, _u: &[f64], sigma: &mut [f64], _grid: &Grid1D<f64>) {
1307 for s in sigma.iter_mut() {
1308 *s = self.sigma;
1309 }
1310 }
1311
1312 fn noise_correlation(&self) -> crate::system::NoiseCorrelation<f64> {
1313 crate::system::NoiseCorrelation::TraceClass {
1314 n_modes: self.n_modes,
1315 decay_rate: self.decay_rate,
1316 }
1317 }
1318 }
1319
1320 #[test]
1321 fn test_spde_trace_class_noise_fixed_step() {
1322 let system = TraceClassHeat {
1323 alpha: 0.01,
1324 sigma: 0.1,
1325 n_modes: 5,
1326 decay_rate: 2.0,
1327 };
1328
1329 let grid = Grid1D::uniform(0.0, 1.0, 22); let u0: Vec<f64> = grid
1331 .interior_points()
1332 .iter()
1333 .map(|&x| (std::f64::consts::PI * x).sin())
1334 .collect();
1335
1336 let options = SpdeOptions::default().dt(0.0001).n_output(20).seed(42);
1337
1338 let result = MolSdeSolver::solve(&system, 0.0, 0.01, &u0, &grid, &options).unwrap();
1339
1340 assert!(result.success);
1341 for y in &result.y {
1342 for &val in y {
1343 assert!(
1344 val.is_finite(),
1345 "Trace-class noise solution should be finite"
1346 );
1347 }
1348 }
1349 }
1350
1351 #[test]
1352 fn test_spde_trace_class_noise_smoother_than_white() {
1353 struct PureNoise;
1357 impl SpdeSystem<f64> for PureNoise {
1358 fn drift(&self, _t: f64, _u: &[f64], du: &mut [f64], _grid: &Grid1D<f64>) {
1359 for d in du.iter_mut() {
1360 *d = 0.0;
1361 }
1362 }
1363 fn diffusion(&self, _t: f64, _u: &[f64], sigma: &mut [f64], _grid: &Grid1D<f64>) {
1364 for s in sigma.iter_mut() {
1365 *s = 1.0;
1366 }
1367 }
1368 }
1369
1370 struct PureTraceNoise;
1371 impl SpdeSystem<f64> for PureTraceNoise {
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 fn noise_correlation(&self) -> crate::system::NoiseCorrelation<f64> {
1383 crate::system::NoiseCorrelation::TraceClass {
1384 n_modes: 3,
1385 decay_rate: 2.0,
1386 }
1387 }
1388 }
1389
1390 let grid = Grid1D::uniform(0.0, 1.0, 32); let n_interior = grid.len() - 2;
1392 let u0 = vec![0.0; n_interior];
1393
1394 let n_traj = 100;
1395 let mut white_tv_total = 0.0;
1396 let mut trace_tv_total = 0.0;
1397
1398 for i in 0..n_traj {
1399 let options = SpdeOptions::default()
1400 .dt(0.001)
1401 .n_output(5)
1402 .seed(2000 + i as u64);
1403
1404 let white_result =
1405 MolSdeSolver::solve(&PureNoise, 0.0, 0.01, &u0, &grid, &options).unwrap();
1406 let trace_result =
1407 MolSdeSolver::solve(&PureTraceNoise, 0.0, 0.01, &u0, &grid, &options).unwrap();
1408
1409 let white_y = white_result.y_final().unwrap();
1410 let trace_y = trace_result.y_final().unwrap();
1411
1412 let white_tv: f64 = white_y.windows(2).map(|w| (w[1] - w[0]).abs()).sum();
1414 let trace_tv: f64 = trace_y.windows(2).map(|w| (w[1] - w[0]).abs()).sum();
1415
1416 white_tv_total += white_tv;
1417 trace_tv_total += trace_tv;
1418 }
1419
1420 assert!(
1422 trace_tv_total < white_tv_total,
1423 "Trace-class noise should be smoother: trace TV={}, white TV={}",
1424 trace_tv_total / n_traj as f64,
1425 white_tv_total / n_traj as f64
1426 );
1427 }
1428}