1use numra_core::Scalar;
28use numra_sde::wiener::{create_wiener, WienerProcess};
29use rand::prelude::*;
30use rand_distr::StandardNormal;
31
32pub trait SpaceTimeNoise<S: Scalar> {
34 fn increment(&mut self, dt: S, n_points: usize, dw: &mut [S]);
41
42 fn reseed(&mut self, seed: u64);
44}
45
46pub struct WhiteNoise<S: Scalar> {
50 wiener: WienerProcess<StdRng>,
51 dx: S,
52 n_wiener: usize,
53}
54
55impl<S: Scalar> WhiteNoise<S> {
56 pub fn new(dx: S, n_points: usize, seed: Option<u64>) -> Self {
63 Self {
64 wiener: create_wiener(n_points, seed),
65 dx,
66 n_wiener: n_points,
67 }
68 }
69}
70
71impl<S: Scalar> SpaceTimeNoise<S> for WhiteNoise<S> {
72 fn increment(&mut self, dt: S, n_points: usize, dw: &mut [S]) {
73 if n_points != self.n_wiener {
75 self.wiener = create_wiener(n_points, None);
76 self.n_wiener = n_points;
77 }
78
79 let scale = S::ONE / self.dx.sqrt();
88
89 let inc = self.wiener.increment::<S>(dt);
90 for i in 0..n_points.min(inc.dw.len()) {
91 dw[i] = inc.dw[i] * scale;
92 }
93 }
94
95 fn reseed(&mut self, seed: u64) {
96 self.wiener = create_wiener(self.n_wiener, Some(seed));
97 }
98}
99
100pub struct ColoredNoise<S: Scalar> {
107 wiener: WienerProcess<StdRng>,
108 dx: S,
109 correlation_length: S,
110 cholesky: Vec<S>,
112 n_cached: usize,
113}
114
115impl<S: Scalar> ColoredNoise<S> {
116 pub fn new(
127 dx: S,
128 correlation_length: S,
129 n_points: usize,
130 seed: Option<u64>,
131 ) -> Result<Self, String> {
132 let cholesky = Self::compute_cholesky(dx, correlation_length, n_points)?;
134
135 Ok(Self {
136 wiener: create_wiener(n_points, seed),
137 dx,
138 correlation_length,
139 cholesky,
140 n_cached: n_points,
141 })
142 }
143
144 fn compute_cholesky(dx: S, correlation_length: S, n: usize) -> Result<Vec<S>, String> {
145 let mut c = vec![S::ZERO; n * n];
147
148 for i in 0..n {
149 for j in 0..n {
150 let dist = dx * S::from_usize(if i > j { i - j } else { j - i });
151 c[i * n + j] = (-dist / correlation_length).exp();
152 }
153 }
154
155 cholesky_decompose(&c, n)
157 }
158}
159
160impl<S: Scalar> SpaceTimeNoise<S> for ColoredNoise<S> {
161 fn increment(&mut self, dt: S, n_points: usize, dw: &mut [S]) {
162 if n_points != self.n_cached {
164 self.cholesky = Self::compute_cholesky(self.dx, self.correlation_length, n_points)
165 .expect("Cholesky decomposition failed on grid resize; correlation matrix is ill-conditioned");
166 self.wiener = create_wiener(n_points, None);
167 self.n_cached = n_points;
168 }
169
170 let inc = self.wiener.increment::<S>(dt);
172
173 for i in 0..n_points {
175 dw[i] = S::ZERO;
176 for j in 0..=i {
177 dw[i] = dw[i] + self.cholesky[i * n_points + j] * inc.dw[j];
178 }
179 }
180 }
181
182 fn reseed(&mut self, seed: u64) {
183 self.wiener = create_wiener(self.n_cached, Some(seed));
184 }
185}
186
187#[allow(dead_code)]
194pub struct TraceClassNoise<S: Scalar> {
195 wiener: WienerProcess<StdRng>,
196 n_modes: usize,
197 eigenvalues: Vec<S>,
198 domain_length: S,
199}
200
201#[allow(dead_code)]
202impl<S: Scalar> TraceClassNoise<S> {
203 pub fn new(n_modes: usize, decay_rate: S, domain_length: S, seed: Option<u64>) -> Self {
211 let eigenvalues: Vec<S> = (1..=n_modes)
213 .map(|k| S::from_usize(k).powf(-decay_rate))
214 .collect();
215
216 Self {
217 wiener: create_wiener(n_modes, seed),
218 n_modes,
219 eigenvalues,
220 domain_length,
221 }
222 }
223}
224
225impl<S: Scalar> SpaceTimeNoise<S> for TraceClassNoise<S> {
226 fn increment(&mut self, dt: S, n_points: usize, dw: &mut [S]) {
227 for i in 0..n_points {
229 dw[i] = S::ZERO;
230 }
231
232 let pi = S::PI;
233
234 let inc = self.wiener.increment::<S>(dt);
236
237 for k in 0..self.n_modes {
239 let mode_idx = S::from_usize(k + 1);
240 let dw_k = inc.dw[k];
241 let sqrt_lambda = self.eigenvalues[k].sqrt();
242
243 let normalization = (S::TWO / self.domain_length).sqrt();
245
246 for i in 0..n_points {
247 let x = self.domain_length * S::from_usize(i) / S::from_usize(n_points - 1);
248 let basis = normalization * (mode_idx * pi * x / self.domain_length).sin();
249 dw[i] = dw[i] + sqrt_lambda * basis * dw_k;
250 }
251 }
252 }
253
254 fn reseed(&mut self, seed: u64) {
255 self.wiener = create_wiener(self.n_modes, Some(seed));
256 }
257}
258
259pub struct ColoredNoiseGenerator<S: Scalar> {
272 correlation_length: S,
273 rng: StdRng,
274}
275
276impl<S: Scalar> ColoredNoiseGenerator<S> {
277 pub fn new(correlation_length: S, seed: u64) -> Self {
283 Self {
284 correlation_length,
285 rng: StdRng::seed_from_u64(seed),
286 }
287 }
288
289 pub fn sample(&mut self, grid: &[S]) -> Result<Vec<S>, String> {
297 let n = grid.len();
298 let l = self.correlation_length;
299
300 let mut c = vec![S::ZERO; n * n];
302 for i in 0..n {
303 for j in 0..n {
304 let dx = grid[i] - grid[j];
305 c[i * n + j] = (-dx * dx / (S::from_f64(2.0) * l * l)).exp();
306 }
307 }
308
309 let chol = cholesky_decompose(&c, n)?;
311
312 let z: Vec<S> = (0..n)
314 .map(|_| S::from_f64(self.rng.sample(StandardNormal)))
315 .collect();
316
317 let mut result = vec![S::ZERO; n];
319 for i in 0..n {
320 for j in 0..=i {
321 result[i] = result[i] + chol[i * n + j] * z[j];
322 }
323 }
324
325 Ok(result)
326 }
327}
328
329pub(crate) fn cholesky_decompose<S: Scalar>(a: &[S], n: usize) -> Result<Vec<S>, String> {
337 let mut l = vec![S::ZERO; n * n];
338 let eps = S::from_f64(1e-12);
339
340 for j in 0..n {
341 let mut sum = S::ZERO;
342 for k in 0..j {
343 sum = sum + l[j * n + k] * l[j * n + k];
344 }
345 let diag = a[j * n + j] - sum;
346 if diag < S::ZERO {
347 return Err(format!(
348 "Correlation matrix not positive definite at index {}: \
349 diagonal element is {}",
350 j,
351 diag.to_f64()
352 ));
353 }
354 let sqrt_diag = diag.sqrt();
355 if sqrt_diag < eps {
356 return Err(format!(
357 "Correlation matrix ill-conditioned at index {}: \
358 diagonal element {} is below threshold. \
359 Try increasing correlation_length.",
360 j,
361 diag.to_f64()
362 ));
363 }
364 l[j * n + j] = sqrt_diag;
365
366 for i in (j + 1)..n {
367 let mut sum = S::ZERO;
368 for k in 0..j {
369 sum = sum + l[i * n + k] * l[j * n + k];
370 }
371 l[i * n + j] = (a[i * n + j] - sum) / l[j * n + j];
372 }
373 }
374
375 Ok(l)
376}
377
378#[cfg(test)]
379mod tests {
380 use super::*;
381
382 #[test]
383 fn test_white_noise_variance() {
384 let dx = 0.1;
385 let dt = 0.01;
386 let n = 10;
387 let n_samples = 1000;
388
389 let mut noise = WhiteNoise::new(dx, n, Some(42));
390 let mut dw = vec![0.0; n];
391 let mut sum_sq = vec![0.0; n];
392
393 for _ in 0..n_samples {
394 noise.increment(dt, n, &mut dw);
395 for i in 0..n {
396 sum_sq[i] += dw[i] * dw[i];
397 }
398 }
399
400 let expected_var = dt / dx;
402 for i in 0..n {
403 let var = sum_sq[i] / n_samples as f64;
404 assert!(
405 (var - expected_var).abs() < 0.05,
406 "Variance at point {}: {}, expected: {}",
407 i,
408 var,
409 expected_var
410 );
411 }
412 }
413
414 #[test]
415 fn test_colored_noise_correlation() {
416 let dx = 0.1;
417 let correlation_length = 0.3;
418 let n = 20;
419 let dt = 0.01;
420 let n_samples = 5000;
421
422 let mut noise = ColoredNoise::new(dx, correlation_length, n, Some(42)).unwrap();
423 let mut dw = vec![0.0; n];
424 let mut cross_sum = 0.0;
425 let mut sum_sq = 0.0;
426
427 for _ in 0..n_samples {
429 noise.increment(dt, n, &mut dw);
430 cross_sum += dw[0] * dw[1];
431 sum_sq += dw[0] * dw[0];
432 }
433
434 let correlation = cross_sum / sum_sq;
435 let expected = (-dx / correlation_length).exp();
436
437 assert!(
439 correlation > 0.0,
440 "Adjacent points should be positively correlated"
441 );
442 assert!(
443 (correlation - expected).abs() < 0.2,
444 "Correlation: {}, expected: {}",
445 correlation,
446 expected
447 );
448 }
449
450 #[test]
451 fn test_trace_class_noise() {
452 let n_modes = 10;
453 let decay_rate = 2.0;
454 let domain_length = 1.0;
455 let n = 50;
456 let dt = 0.01;
457
458 let mut noise = TraceClassNoise::new(n_modes, decay_rate, domain_length, Some(42));
459 let mut dw = vec![0.0; n];
460
461 noise.increment(dt, n, &mut dw);
462
463 let mut total_variation = 0.0;
466 for i in 1..n {
467 total_variation += (dw[i] - dw[i - 1]).abs();
468 }
469
470 assert!(total_variation.is_finite());
472 }
473
474 #[test]
475 fn test_colored_noise_generation() {
476 use crate::noise::ColoredNoiseGenerator;
477
478 let correlation_length = 0.1_f64;
479 let grid: Vec<f64> = (0..20).map(|i| i as f64 / 19.0).collect();
480 let n = grid.len();
481
482 let mut generator = ColoredNoiseGenerator::new(correlation_length, 42);
483
484 let n_samples = 2000;
486 let mut samples = Vec::with_capacity(n_samples);
487 for _ in 0..n_samples {
488 samples.push(generator.sample(&grid).unwrap());
489 }
490
491 let mean_0: f64 = samples.iter().map(|s| s[0]).sum::<f64>() / n_samples as f64;
496 let mean_1: f64 = samples.iter().map(|s| s[1]).sum::<f64>() / n_samples as f64;
497 let var_0: f64 =
498 samples.iter().map(|s| (s[0] - mean_0).powi(2)).sum::<f64>() / n_samples as f64;
499 let var_1: f64 =
500 samples.iter().map(|s| (s[1] - mean_1).powi(2)).sum::<f64>() / n_samples as f64;
501 let cov_01: f64 = samples
502 .iter()
503 .map(|s| (s[0] - mean_0) * (s[1] - mean_1))
504 .sum::<f64>()
505 / n_samples as f64;
506 let corr_01 = cov_01 / (var_0.sqrt() * var_1.sqrt());
507
508 assert!(corr_01 > 0.5, "Adjacent correlation too low: {}", corr_01);
510
511 let mean_far: f64 = samples.iter().map(|s| s[n - 1]).sum::<f64>() / n_samples as f64;
513 let var_far: f64 = samples
514 .iter()
515 .map(|s| (s[n - 1] - mean_far).powi(2))
516 .sum::<f64>()
517 / n_samples as f64;
518 let cov_far: f64 = samples
519 .iter()
520 .map(|s| (s[0] - mean_0) * (s[n - 1] - mean_far))
521 .sum::<f64>()
522 / n_samples as f64;
523 let corr_far = cov_far / (var_0.sqrt() * var_far.sqrt());
524
525 assert!(
526 corr_01 > corr_far.abs(),
527 "Nearby correlation should be higher than far correlation"
528 );
529 }
530
531 #[test]
532 fn test_cholesky_rejects_non_positive_definite() {
533 let n = 2;
535 let matrix = vec![1.0_f64, 2.0, 2.0, 1.0];
536 let result = cholesky_decompose::<f64>(&matrix, n);
537 assert!(
538 result.is_err(),
539 "Should reject non-positive-definite matrix"
540 );
541 let msg = result.unwrap_err();
542 assert!(
543 msg.contains("not positive definite"),
544 "Error message: {}",
545 msg
546 );
547 }
548
549 #[test]
550 fn test_cholesky_rejects_ill_conditioned() {
551 let n = 3;
553 let mut matrix = vec![0.0_f64; n * n];
554 matrix[0] = 1.0;
555 matrix[4] = 1e-30; matrix[8] = 1.0;
557 let result = cholesky_decompose::<f64>(&matrix, n);
558 assert!(result.is_err(), "Should reject ill-conditioned matrix");
559 let msg = result.unwrap_err();
560 assert!(msg.contains("ill-conditioned"), "Error message: {}", msg);
561 }
562
563 #[test]
564 fn test_cholesky_accepts_valid_matrix() {
565 let n = 3;
567 let mut matrix = vec![0.0_f64; n * n];
568 matrix[0] = 1.0;
569 matrix[4] = 1.0;
570 matrix[8] = 1.0;
571 let result = cholesky_decompose::<f64>(&matrix, n);
572 assert!(result.is_ok(), "Should accept identity matrix");
573 let l = result.unwrap();
574 assert!((l[0] - 1.0).abs() < 1e-10);
576 assert!((l[4] - 1.0).abs() < 1e-10);
577 assert!((l[8] - 1.0).abs() < 1e-10);
578 }
579
580 #[test]
581 fn test_colored_noise_new_returns_ok_for_valid_params() {
582 let result = ColoredNoise::<f64>::new(0.1, 0.3, 10, Some(42));
583 assert!(
584 result.is_ok(),
585 "Should accept reasonable correlation parameters"
586 );
587 }
588
589 #[test]
590 fn test_cholesky_decompose_recovers_factor() {
591 let n = 3;
593 let a = vec![4.0_f64, 2.0, 1.0, 2.0, 5.0, 3.0, 1.0, 3.0, 6.0];
594 let l = cholesky_decompose::<f64>(&a, n).unwrap();
595
596 for i in 0..n {
598 for j in 0..n {
599 let mut sum = 0.0;
600 for k in 0..n {
601 sum += l[i * n + k] * l[j * n + k];
602 }
603 assert!(
604 (sum - a[i * n + j]).abs() < 1e-10,
605 "Reconstruction error at ({},{}): {} vs {}",
606 i,
607 j,
608 sum,
609 a[i * n + j]
610 );
611 }
612 }
613 }
614}