scirs2_integrate/
verification.rs

1//! Method of Manufactured Solutions (MMS) toolkit for verification
2//!
3//! This module provides tools for verifying numerical methods using the Method of
4//! Manufactured Solutions. MMS is a powerful technique for code verification where
5//! an exact solution is chosen, the governing equation is modified by adding source
6//! terms, and the numerical solution is compared against the exact solution.
7//!
8//! # Examples
9//!
10//! ## ODE Verification
11//! ```
12//! use scirs2_integrate::verification::{MMSODEProblem, polynomial_solution};
13//! use scirs2_core::ndarray::Array1;
14//!
15//! // Create a manufactured ODE problem with polynomial solution
16//! let exact_solution = polynomial_solution(vec![1.0, 2.0, 3.0]); // 1 + 2t + 3t²
17//! let problem = MMSODEProblem::new(exact_solution, [0.0, 1.0]);
18//!
19//! // The manufactured source term is automatically computed
20//! // Solve numerically and verify order of accuracy
21//! ```
22//!
23//! ## PDE Verification  
24//! ```
25//! use scirs2_integrate::verification::{MMSPDEProblem, trigonometric_solution_2d};
26//!
27//! // Create manufactured 2D Poisson problem
28//! let exact = trigonometric_solution_2d(1.0, 2.0); // sin(x) * cos(2y)
29//! let problem = MMSPDEProblem::new_poisson_2d(exact, [0.0, 1.0], [0.0, 1.0]);
30//! ```
31
32use crate::common::IntegrateFloat;
33use crate::error::{IntegrateError, IntegrateResult};
34use scirs2_core::ndarray::{ArrayView1, ArrayView2};
35use std::f64::consts::PI;
36use std::fmt;
37
38/// Trait for exact solutions in MMS problems
39pub trait ExactSolution<F: IntegrateFloat> {
40    /// Evaluate the exact solution at given point(s)
41    fn evaluate(&self, coordinates: &[F]) -> F;
42
43    /// Evaluate the first derivative with respect to specified variable
44    fn derivative(&self, coordinates: &[F], variable: usize) -> F;
45
46    /// Evaluate the second derivative with respect to specified variable
47    fn second_derivative(&self, coordinates: &[F], variable: usize) -> F;
48
49    /// Evaluate mixed partial derivatives (for PDEs)
50    fn mixed_derivative(&self, coordinates: &[F], var1: usize, var2: usize) -> F {
51        // Default implementation using finite differences
52        let h = F::from(1e-8).unwrap();
53        let mut coords_plus = coordinates.to_vec();
54        let mut coords_minus = coordinates.to_vec();
55
56        coords_plus[var2] += h;
57        coords_minus[var2] -= h;
58
59        let deriv_plus = self.derivative(&coords_plus, var1);
60        let deriv_minus = self.derivative(&coords_minus, var1);
61
62        (deriv_plus - deriv_minus) / (F::from(2.0).unwrap() * h)
63    }
64
65    /// Get the dimensionality of the problem
66    fn dimension(&self) -> usize;
67}
68
69/// Polynomial exact solution for ODE problems
70#[derive(Debug, Clone)]
71pub struct PolynomialSolution<F: IntegrateFloat> {
72    coefficients: Vec<F>,
73}
74
75impl<F: IntegrateFloat> PolynomialSolution<F> {
76    /// Create a new polynomial solution: sum(coeff[i] * t^i)
77    pub fn new(coefficients: Vec<F>) -> Self {
78        Self { coefficients }
79    }
80}
81
82impl<F: IntegrateFloat> ExactSolution<F> for PolynomialSolution<F> {
83    fn evaluate(&self, coordinates: &[F]) -> F {
84        let t = coordinates[0];
85        let mut result = F::zero();
86        let mut t_power = F::one();
87
88        for &coeff in &self.coefficients {
89            result += coeff * t_power;
90            t_power *= t;
91        }
92
93        result
94    }
95
96    fn derivative(&self, coordinates: &[F], variable: usize) -> F {
97        let t = coordinates[0];
98        let mut result = F::zero();
99        let mut t_power = F::one();
100
101        for (i, &coeff) in self.coefficients.iter().enumerate().skip(1) {
102            result += F::from(i).unwrap() * coeff * t_power;
103            t_power *= t;
104        }
105
106        result
107    }
108
109    fn second_derivative(&self, coordinates: &[F], variable: usize) -> F {
110        let t = coordinates[0];
111        let mut result = F::zero();
112        let mut t_power = F::one();
113
114        for (i, &coeff) in self.coefficients.iter().enumerate().skip(2) {
115            let factor = F::from(i * (i - 1)).unwrap();
116            result += factor * coeff * t_power;
117            t_power *= t;
118        }
119
120        result
121    }
122
123    fn dimension(&self) -> usize {
124        1
125    }
126}
127
128/// Trigonometric exact solution for 2D PDE problems  
129#[derive(Debug, Clone)]
130pub struct TrigonometricSolution2D<F: IntegrateFloat> {
131    freq_x: F,
132    freq_y: F,
133    phase_x: F,
134    phase_y: F,
135}
136
137impl<F: IntegrateFloat> TrigonometricSolution2D<F> {
138    /// Create sin(freq_x * x + phase_x) * cos(freq_y * y + phase_y)
139    pub fn new(freq_x: F, freq_y: F, phase_x: F, phase_y: F) -> Self {
140        Self {
141            freq_x,
142            freq_y,
143            phase_x,
144            phase_y,
145        }
146    }
147
148    /// Create sin(freq_x * x) * cos(freq_y * y)
149    pub fn simple(freq_x: F, freqy: F) -> Self {
150        Self::new(freq_x, freqy, F::zero(), F::zero())
151    }
152}
153
154impl<F: IntegrateFloat> ExactSolution<F> for TrigonometricSolution2D<F> {
155    fn evaluate(&self, coordinates: &[F]) -> F {
156        let x = coordinates[0];
157        let y = coordinates[1];
158        (self.freq_x * x + self.phase_x).sin() * (self.freq_y * y + self.phase_y).cos()
159    }
160
161    fn derivative(&self, coordinates: &[F], variable: usize) -> F {
162        let x = coordinates[0];
163        let y = coordinates[1];
164
165        match variable {
166            0 => {
167                self.freq_x
168                    * (self.freq_x * x + self.phase_x).cos()
169                    * (self.freq_y * y + self.phase_y).cos()
170            }
171            1 => {
172                -self.freq_y
173                    * (self.freq_x * x + self.phase_x).sin()
174                    * (self.freq_y * y + self.phase_y).sin()
175            }
176            _ => F::zero(),
177        }
178    }
179
180    fn second_derivative(&self, coordinates: &[F], variable: usize) -> F {
181        let x = coordinates[0];
182        let y = coordinates[1];
183
184        match variable {
185            0 => {
186                -self.freq_x
187                    * self.freq_x
188                    * (self.freq_x * x + self.phase_x).sin()
189                    * (self.freq_y * y + self.phase_y).cos()
190            }
191            1 => {
192                -self.freq_y
193                    * self.freq_y
194                    * (self.freq_x * x + self.phase_x).sin()
195                    * (self.freq_y * y + self.phase_y).cos()
196            }
197            _ => F::zero(),
198        }
199    }
200
201    fn dimension(&self) -> usize {
202        2
203    }
204}
205
206/// Manufactured ODE problem for verification
207#[derive(Debug)]
208pub struct MMSODEProblem<F: IntegrateFloat, S: ExactSolution<F>> {
209    exact_solution: S,
210    time_span: [F; 2],
211    _phantom: std::marker::PhantomData<F>,
212}
213
214impl<F: IntegrateFloat, S: ExactSolution<F>> MMSODEProblem<F, S> {
215    /// Create a new manufactured ODE problem
216    pub fn new(exact_solution: S, timespan: [F; 2]) -> Self {
217        Self {
218            exact_solution,
219            time_span: timespan,
220            _phantom: std::marker::PhantomData,
221        }
222    }
223
224    /// Get the manufactured source term for: y' = f(t, y) + source(t)
225    /// where f is the original RHS and source is computed from exact solution
226    pub fn source_term(&self, t: F) -> F {
227        let coords = [t];
228        // For y' = f(t, y), source = y_exact'(t) - f(t, y_exact(t))
229        // If f is linear: y' = a*y + b, then source = y_exact'(t) - a*y_exact(t) - b
230        // For autonomous f: y' = f(y), source = y_exact'(t) - f(y_exact(t))
231        self.exact_solution.derivative(&coords, 0)
232    }
233
234    /// Get initial condition from exact solution
235    pub fn initial_condition(&self) -> F {
236        self.exact_solution.evaluate(&[self.time_span[0]])
237    }
238
239    /// Evaluate exact solution at given time
240    pub fn exact_at(&self, t: F) -> F {
241        self.exact_solution.evaluate(&[t])
242    }
243
244    /// Get time span
245    pub fn time_span(&self) -> [F; 2] {
246        self.time_span
247    }
248}
249
250/// Manufactured PDE problem for verification  
251#[derive(Debug)]
252pub struct MMSPDEProblem<F: IntegrateFloat, S: ExactSolution<F>> {
253    exact_solution: S,
254    domain_x: [F; 2],
255    domain_y: [F; 2],
256    domain_z: Option<[F; 2]>,
257    pde_type: PDEType,
258    parameters: PDEParameters<F>,
259    _phantom: std::marker::PhantomData<F>,
260}
261
262/// Parameters for different PDE types
263#[derive(Debug, Clone)]
264pub struct PDEParameters<F: IntegrateFloat> {
265    /// Diffusion coefficient (for diffusion/heat equation)
266    pub diffusion_coeff: F,
267    /// Wave speed (for wave equation)
268    pub wave_speed: F,
269    /// Advection velocity (for advection-diffusion)
270    pub advection_velocity: Vec<F>,
271    /// Reaction coefficient (for reaction-diffusion)
272    pub reaction_coeff: F,
273    /// Helmholtz parameter (for Helmholtz equation)
274    pub helmholtz_k: F,
275}
276
277impl<F: IntegrateFloat> Default for PDEParameters<F> {
278    fn default() -> Self {
279        Self {
280            diffusion_coeff: F::one(),
281            wave_speed: F::one(),
282            advection_velocity: vec![F::zero()],
283            reaction_coeff: F::zero(),
284            helmholtz_k: F::one(),
285        }
286    }
287}
288
289#[derive(Debug, Clone)]
290pub enum PDEType {
291    Poisson2D,
292    Poisson3D,
293    Diffusion2D,
294    Diffusion3D,
295    Wave2D,
296    Wave3D,
297    AdvectionDiffusion2D,
298    AdvectionDiffusion3D,
299    Helmholtz2D,
300    Helmholtz3D,
301}
302
303impl<F: IntegrateFloat, S: ExactSolution<F>> MMSPDEProblem<F, S> {
304    /// Create manufactured 2D Poisson problem: -∇²u = f
305    pub fn new_poisson_2d(exact_solution: S, domain_x: [F; 2], domainy: [F; 2]) -> Self {
306        Self {
307            exact_solution,
308            domain_x,
309            domain_y: domainy,
310            domain_z: None,
311            pde_type: PDEType::Poisson2D,
312            parameters: PDEParameters::default(),
313            _phantom: std::marker::PhantomData,
314        }
315    }
316
317    /// Create manufactured 3D Poisson problem: -∇²u = f
318    pub fn new_poisson_3d(
319        exact_solution: S,
320        domain_x: [F; 2],
321        domain_y: [F; 2],
322        domain_z: [F; 2],
323    ) -> Self {
324        Self {
325            exact_solution,
326            domain_x,
327            domain_y,
328            domain_z: Some(domain_z),
329            pde_type: PDEType::Poisson3D,
330            parameters: PDEParameters::default(),
331            _phantom: std::marker::PhantomData,
332        }
333    }
334
335    /// Create manufactured 2D diffusion problem: ∂u/∂t = α∇²u + f
336    pub fn new_diffusion_2d(
337        exact_solution: S,
338        domain_x: [F; 2],
339        domain_y: [F; 2],
340        diffusion_coeff: F,
341    ) -> Self {
342        let mut params = PDEParameters::default();
343        params.diffusion_coeff = diffusion_coeff;
344        Self {
345            exact_solution,
346            domain_x,
347            domain_y,
348            domain_z: None,
349            pde_type: PDEType::Diffusion2D,
350            parameters: params,
351            _phantom: std::marker::PhantomData,
352        }
353    }
354
355    /// Create manufactured 2D wave problem: ∂²u/∂t² = c²∇²u + f
356    pub fn new_wave_2d(
357        exact_solution: S,
358        domain_x: [F; 2],
359        domain_y: [F; 2],
360        wave_speed: F,
361    ) -> Self {
362        let mut params = PDEParameters::default();
363        params.wave_speed = wave_speed;
364        Self {
365            exact_solution,
366            domain_x,
367            domain_y,
368            domain_z: None,
369            pde_type: PDEType::Wave2D,
370            parameters: params,
371            _phantom: std::marker::PhantomData,
372        }
373    }
374
375    /// Create manufactured 2D Helmholtz problem: ∇²u + k²u = f
376    pub fn new_helmholtz_2d(exact_solution: S, domain_x: [F; 2], domainy: [F; 2], k: F) -> Self {
377        let mut params = PDEParameters::default();
378        params.helmholtz_k = k;
379        Self {
380            exact_solution,
381            domain_x,
382            domain_y: domainy,
383            domain_z: None,
384            pde_type: PDEType::Helmholtz2D,
385            parameters: params,
386            _phantom: std::marker::PhantomData,
387        }
388    }
389
390    /// Get manufactured source term for the PDE
391    pub fn source_term(&self, coordinates: &[F]) -> F {
392        match self.pde_type {
393            PDEType::Poisson2D => {
394                // For -∇²u = f, source f = -∇²u_exact
395                -(self.exact_solution.second_derivative(coordinates, 0)
396                    + self.exact_solution.second_derivative(coordinates, 1))
397            }
398            PDEType::Poisson3D => {
399                // For -∇²u = f in 3D
400                -(self.exact_solution.second_derivative(coordinates, 0)
401                    + self.exact_solution.second_derivative(coordinates, 1)
402                    + self.exact_solution.second_derivative(coordinates, 2))
403            }
404            PDEType::Diffusion2D => {
405                // For ∂u/∂t = α∇²u + f
406                // source f = ∂u_exact/∂t - α∇²u_exact
407                let time_dim = coordinates.len() - 1;
408                let spatial_laplacian = self.exact_solution.second_derivative(coordinates, 0)
409                    + self.exact_solution.second_derivative(coordinates, 1);
410                self.exact_solution.derivative(coordinates, time_dim)
411                    - self.parameters.diffusion_coeff * spatial_laplacian
412            }
413            PDEType::Wave2D => {
414                // For ∂²u/∂t² = c²∇²u + f
415                // source f = ∂²u_exact/∂t² - c²∇²u_exact
416                // Simplified: assume time coordinate is last
417                let spatial_laplacian = self.exact_solution.second_derivative(coordinates, 0)
418                    + self.exact_solution.second_derivative(coordinates, 1);
419                // Second time derivative would need higher-order implementation
420                F::zero()
421                    - self.parameters.wave_speed * self.parameters.wave_speed * spatial_laplacian
422            }
423            PDEType::Helmholtz2D => {
424                // For ∇²u + k²u = f
425                // source f = ∇²u_exact + k²u_exact
426                let laplacian = self.exact_solution.second_derivative(coordinates, 0)
427                    + self.exact_solution.second_derivative(coordinates, 1);
428                laplacian
429                    + self.parameters.helmholtz_k
430                        * self.parameters.helmholtz_k
431                        * self.exact_solution.evaluate(coordinates)
432            }
433            PDEType::AdvectionDiffusion2D => {
434                // For ∂u/∂t + v·∇u = α∇²u + f
435                // source f = ∂u_exact/∂t + v·∇u_exact - α∇²u_exact
436                let time_dim = coordinates.len() - 1;
437                let time_deriv = self.exact_solution.derivative(coordinates, time_dim);
438                let spatial_laplacian = self.exact_solution.second_derivative(coordinates, 0)
439                    + self.exact_solution.second_derivative(coordinates, 1);
440
441                let mut advection_term = F::zero();
442                for (i, &v_i) in self
443                    .parameters
444                    .advection_velocity
445                    .iter()
446                    .enumerate()
447                    .take(2)
448                {
449                    advection_term += v_i * self.exact_solution.derivative(coordinates, i);
450                }
451
452                time_deriv + advection_term - self.parameters.diffusion_coeff * spatial_laplacian
453            }
454            _ => F::zero(),
455        }
456    }
457
458    /// Get boundary condition from exact solution (2D)
459    pub fn boundary_condition(&self, x: F, y: F) -> F {
460        self.exact_solution.evaluate(&[x, y])
461    }
462
463    /// Get boundary condition from exact solution (3D)
464    pub fn boundary_condition_3d(&self, x: F, y: F, z: F) -> F {
465        self.exact_solution.evaluate(&[x, y, z])
466    }
467
468    /// Evaluate exact solution (2D)
469    pub fn exact_at(&self, x: F, y: F) -> F {
470        self.exact_solution.evaluate(&[x, y])
471    }
472
473    /// Evaluate exact solution (3D)
474    pub fn exact_at_3d(&self, x: F, y: F, z: F) -> F {
475        self.exact_solution.evaluate(&[x, y, z])
476    }
477
478    /// Get domain bounds (2D)
479    pub fn domain(&self) -> ([F; 2], [F; 2]) {
480        (self.domain_x, self.domain_y)
481    }
482
483    /// Get domain bounds (3D)
484    pub fn domain_3d(&self) -> ([F; 2], [F; 2], [F; 2]) {
485        (
486            self.domain_x,
487            self.domain_y,
488            self.domain_z.unwrap_or([F::zero(), F::one()]),
489        )
490    }
491
492    /// Get PDE parameters
493    pub fn parameters(&self) -> &PDEParameters<F> {
494        &self.parameters
495    }
496
497    /// Check if problem is 3D
498    pub fn is_3d(&self) -> bool {
499        self.domain_z.is_some()
500    }
501}
502
503/// Convergence analysis results
504#[derive(Debug, Clone)]
505pub struct ConvergenceAnalysis<F: IntegrateFloat> {
506    /// Grid sizes or time steps used
507    pub grid_sizes: Vec<F>,
508    /// Error norms for each grid size
509    pub errors: Vec<F>,
510    /// Estimated order of accuracy
511    pub order: F,
512    /// Confidence interval for the order
513    pub order_confidence: (F, F),
514}
515
516impl<F: IntegrateFloat> ConvergenceAnalysis<F> {
517    /// Compute order of accuracy from grid sizes and errors
518    pub fn compute_order(_gridsizes: Vec<F>, errors: Vec<F>) -> IntegrateResult<Self> {
519        if _gridsizes.len() != errors.len() || _gridsizes.len() < 2 {
520            return Err(IntegrateError::ValueError(
521                "Need at least 2 points for convergence analysis".to_string(),
522            ));
523        }
524
525        // Use least squares to fit log(error) = log(C) + p*log(h)
526        let n = _gridsizes.len();
527        let mut sum_log_h = F::zero();
528        let mut sum_log_e = F::zero();
529        let mut sum_log_h_sq = F::zero();
530        let mut sum_log_h_log_e = F::zero();
531
532        for (h, e) in _gridsizes.iter().zip(errors.iter()) {
533            if *e <= F::zero() || *h <= F::zero() {
534                return Err(IntegrateError::ValueError(
535                    "Grid sizes and errors must be positive".to_string(),
536                ));
537            }
538
539            let log_h = h.ln();
540            let log_e = e.ln();
541
542            sum_log_h += log_h;
543            sum_log_e += log_e;
544            sum_log_h_sq += log_h * log_h;
545            sum_log_h_log_e += log_h * log_e;
546        }
547
548        let n_f = F::from(n).unwrap();
549        let denominator = n_f * sum_log_h_sq - sum_log_h * sum_log_h;
550
551        if denominator.abs() < F::from(1e-12).unwrap() {
552            return Err(IntegrateError::ComputationError(
553                "Cannot compute order - insufficient variation in grid sizes".to_string(),
554            ));
555        }
556
557        let order = (n_f * sum_log_h_log_e - sum_log_h * sum_log_e) / denominator;
558
559        // Simple confidence interval (±0.1 for demonstration)
560        let confidence_delta = F::from(0.1).unwrap();
561        let order_confidence = (order - confidence_delta, order + confidence_delta);
562
563        Ok(Self {
564            grid_sizes: _gridsizes,
565            errors,
566            order,
567            order_confidence,
568        })
569    }
570
571    /// Check if the computed order matches expected order within tolerance
572    pub fn verify_order(&self, expectedorder: F, tolerance: F) -> bool {
573        (self.order - expectedorder).abs() <= tolerance
574    }
575}
576
577impl<F: IntegrateFloat> fmt::Display for ConvergenceAnalysis<F> {
578    fn fmt(&self, f: &mut fmt::Formatter<'_>) -> fmt::Result {
579        writeln!(f, "Convergence Analysis Results:")?;
580        writeln!(f, "Grid Size    Error")?;
581        writeln!(f, "─────────────────────")?;
582
583        for (h, e) in self.grid_sizes.iter().zip(self.errors.iter()) {
584            writeln!(f, "{h:9.2e}  {e:9.2e}")?;
585        }
586
587        writeln!(f, "─────────────────────")?;
588        writeln!(f, "Estimated order: {:.3}", self.order)?;
589        writeln!(
590            f,
591            "95% confidence: ({:.3}, {:.3})",
592            self.order_confidence.0, self.order_confidence.1
593        )?;
594
595        Ok(())
596    }
597}
598
599/// Error analysis utilities
600pub struct ErrorAnalysis;
601
602impl ErrorAnalysis {
603    /// Compute L2 norm of error
604    pub fn l2_norm<F: IntegrateFloat>(
605        exact: ArrayView1<F>,
606        numerical: ArrayView1<F>,
607    ) -> IntegrateResult<F> {
608        if exact.len() != numerical.len() {
609            return Err(IntegrateError::ValueError(
610                "Arrays must have same length".to_string(),
611            ));
612        }
613
614        let mut sum_sq = F::zero();
615        for (e, n) in exact.iter().zip(numerical.iter()) {
616            let diff = *e - *n;
617            sum_sq += diff * diff;
618        }
619
620        Ok((sum_sq / F::from(exact.len()).unwrap()).sqrt())
621    }
622
623    /// Compute maximum norm of error
624    pub fn max_norm<F: IntegrateFloat>(
625        exact: ArrayView1<F>,
626        numerical: ArrayView1<F>,
627    ) -> IntegrateResult<F> {
628        if exact.len() != numerical.len() {
629            return Err(IntegrateError::ValueError(
630                "Arrays must have same length".to_string(),
631            ));
632        }
633
634        let mut max_error = F::zero();
635        for (e, n) in exact.iter().zip(numerical.iter()) {
636            let diff = (*e - *n).abs();
637            if diff > max_error {
638                max_error = diff;
639            }
640        }
641
642        Ok(max_error)
643    }
644
645    /// Compute L2 norm of error for 2D arrays
646    pub fn l2_norm_2d<F: IntegrateFloat>(
647        exact: ArrayView2<F>,
648        numerical: ArrayView2<F>,
649    ) -> IntegrateResult<F> {
650        if exact.shape() != numerical.shape() {
651            return Err(IntegrateError::ValueError(
652                "Arrays must have same shape".to_string(),
653            ));
654        }
655
656        let mut sum_sq = F::zero();
657        let mut count = 0;
658
659        for (e, n) in exact.iter().zip(numerical.iter()) {
660            let diff = *e - *n;
661            sum_sq += diff * diff;
662            count += 1;
663        }
664
665        Ok((sum_sq / F::from(count).unwrap()).sqrt())
666    }
667}
668
669/// Convenience functions for creating common exact solutions
670#[allow(dead_code)]
671pub fn polynomial_solution<F: IntegrateFloat>(coefficients: Vec<F>) -> PolynomialSolution<F> {
672    PolynomialSolution::new(coefficients)
673}
674
675#[allow(dead_code)]
676pub fn trigonometric_solution_2d<F: IntegrateFloat>(
677    freq_x: F,
678    freq_y: F,
679) -> TrigonometricSolution2D<F> {
680    TrigonometricSolution2D::simple(freq_x, freq_y)
681}
682
683/// Exponential exact solution for problems with exponential behavior
684#[derive(Debug, Clone)]
685pub struct ExponentialSolution<F: IntegrateFloat> {
686    amplitude: F,
687    decay_rate: F,
688    phase: F,
689}
690
691impl<F: IntegrateFloat> ExponentialSolution<F> {
692    /// Create exp(decay_rate * t + phase) solution
693    pub fn new(_amplitude: F, decayrate: F, phase: F) -> Self {
694        Self {
695            amplitude: _amplitude,
696            decay_rate: decayrate,
697            phase,
698        }
699    }
700
701    /// Create simple exp(decay_rate * t) solution
702    pub fn simple(amplitude: F, decayrate: F) -> Self {
703        Self::new(amplitude, decayrate, F::zero())
704    }
705}
706
707impl<F: IntegrateFloat> ExactSolution<F> for ExponentialSolution<F> {
708    fn evaluate(&self, coordinates: &[F]) -> F {
709        let t = coordinates[0];
710        self.amplitude * (self.decay_rate * t + self.phase).exp()
711    }
712
713    fn derivative(&self, coordinates: &[F], variable: usize) -> F {
714        let t = coordinates[0];
715        self.amplitude * self.decay_rate * (self.decay_rate * t + self.phase).exp()
716    }
717
718    fn second_derivative(&self, coordinates: &[F], variable: usize) -> F {
719        let t = coordinates[0];
720        self.amplitude
721            * self.decay_rate
722            * self.decay_rate
723            * (self.decay_rate * t + self.phase).exp()
724    }
725
726    fn dimension(&self) -> usize {
727        1
728    }
729}
730
731/// Combined solution: polynomial + trigonometric + exponential
732#[derive(Debug, Clone)]
733pub struct CombinedSolution<F: IntegrateFloat> {
734    polynomial: Option<PolynomialSolution<F>>,
735    trigonometric: Option<TrigonometricSolution2D<F>>,
736    exponential: Option<ExponentialSolution<F>>,
737    dimension: usize,
738}
739
740impl<F: IntegrateFloat> CombinedSolution<F> {
741    /// Create a new combined solution
742    pub fn new(dimension: usize) -> Self {
743        Self {
744            polynomial: None,
745            trigonometric: None,
746            exponential: None,
747            dimension,
748        }
749    }
750
751    /// Add polynomial component
752    pub fn with_polynomial(mut self, poly: PolynomialSolution<F>) -> Self {
753        self.polynomial = Some(poly);
754        self
755    }
756
757    /// Add trigonometric component (for 2D problems)
758    pub fn with_trigonometric(mut self, trig: TrigonometricSolution2D<F>) -> Self {
759        self.trigonometric = Some(trig);
760        self
761    }
762
763    /// Add exponential component
764    pub fn with_exponential(mut self, exp: ExponentialSolution<F>) -> Self {
765        self.exponential = Some(exp);
766        self
767    }
768}
769
770impl<F: IntegrateFloat> ExactSolution<F> for CombinedSolution<F> {
771    fn evaluate(&self, coordinates: &[F]) -> F {
772        let mut result = F::zero();
773
774        if let Some(ref poly) = self.polynomial {
775            result += poly.evaluate(coordinates);
776        }
777
778        if let Some(ref trig) = self.trigonometric {
779            if coordinates.len() >= 2 {
780                result += trig.evaluate(coordinates);
781            }
782        }
783
784        if let Some(ref exp) = self.exponential {
785            result += exp.evaluate(coordinates);
786        }
787
788        result
789    }
790
791    fn derivative(&self, coordinates: &[F], variable: usize) -> F {
792        let mut result = F::zero();
793
794        if let Some(ref poly) = self.polynomial {
795            result += poly.derivative(coordinates, variable);
796        }
797
798        if let Some(ref trig) = self.trigonometric {
799            if coordinates.len() >= 2 {
800                result += trig.derivative(coordinates, variable);
801            }
802        }
803
804        if let Some(ref exp) = self.exponential {
805            result += exp.derivative(coordinates, variable);
806        }
807
808        result
809    }
810
811    fn second_derivative(&self, coordinates: &[F], variable: usize) -> F {
812        let mut result = F::zero();
813
814        if let Some(ref poly) = self.polynomial {
815            result += poly.second_derivative(coordinates, variable);
816        }
817
818        if let Some(ref trig) = self.trigonometric {
819            if coordinates.len() >= 2 {
820                result += trig.second_derivative(coordinates, variable);
821            }
822        }
823
824        if let Some(ref exp) = self.exponential {
825            result += exp.second_derivative(coordinates, variable);
826        }
827
828        result
829    }
830
831    fn dimension(&self) -> usize {
832        self.dimension
833    }
834}
835
836/// 3D trigonometric solution for 3D PDE problems
837#[derive(Debug, Clone)]
838pub struct TrigonometricSolution3D<F: IntegrateFloat> {
839    freq_x: F,
840    freq_y: F,
841    freq_z: F,
842    phase_x: F,
843    phase_y: F,
844    phase_z: F,
845}
846
847impl<F: IntegrateFloat> TrigonometricSolution3D<F> {
848    /// Create sin(freq_x * x + phase_x) * cos(freq_y * y + phase_y) * sin(freq_z * z + phase_z)
849    #[allow(clippy::too_many_arguments)]
850    pub fn new(freq_x: F, freq_y: F, freq_z: F, phase_x: F, phase_y: F, phase_z: F) -> Self {
851        Self {
852            freq_x,
853            freq_y,
854            freq_z,
855            phase_x,
856            phase_y,
857            phase_z,
858        }
859    }
860
861    /// Create sin(freq_x * x) * cos(freq_y * y) * sin(freq_z * z)
862    pub fn simple(freq_x: F, freq_y: F, freq_z: F) -> Self {
863        Self::new(freq_x, freq_y, freq_z, F::zero(), F::zero(), F::zero())
864    }
865}
866
867impl<F: IntegrateFloat> ExactSolution<F> for TrigonometricSolution3D<F> {
868    fn evaluate(&self, coordinates: &[F]) -> F {
869        let x = coordinates[0];
870        let y = coordinates[1];
871        let z = coordinates[2];
872        (self.freq_x * x + self.phase_x).sin()
873            * (self.freq_y * y + self.phase_y).cos()
874            * (self.freq_z * z + self.phase_z).sin()
875    }
876
877    fn derivative(&self, coordinates: &[F], variable: usize) -> F {
878        let x = coordinates[0];
879        let y = coordinates[1];
880        let z = coordinates[2];
881
882        match variable {
883            0 => {
884                self.freq_x
885                    * (self.freq_x * x + self.phase_x).cos()
886                    * (self.freq_y * y + self.phase_y).cos()
887                    * (self.freq_z * z + self.phase_z).sin()
888            }
889            1 => {
890                -self.freq_y
891                    * (self.freq_x * x + self.phase_x).sin()
892                    * (self.freq_y * y + self.phase_y).sin()
893                    * (self.freq_z * z + self.phase_z).sin()
894            }
895            2 => {
896                self.freq_z
897                    * (self.freq_x * x + self.phase_x).sin()
898                    * (self.freq_y * y + self.phase_y).cos()
899                    * (self.freq_z * z + self.phase_z).cos()
900            }
901            _ => F::zero(),
902        }
903    }
904
905    fn second_derivative(&self, coordinates: &[F], variable: usize) -> F {
906        let x = coordinates[0];
907        let y = coordinates[1];
908        let z = coordinates[2];
909
910        match variable {
911            0 => {
912                -self.freq_x
913                    * self.freq_x
914                    * (self.freq_x * x + self.phase_x).sin()
915                    * (self.freq_y * y + self.phase_y).cos()
916                    * (self.freq_z * z + self.phase_z).sin()
917            }
918            1 => {
919                -self.freq_y
920                    * self.freq_y
921                    * (self.freq_x * x + self.phase_x).sin()
922                    * (self.freq_y * y + self.phase_y).cos()
923                    * (self.freq_z * z + self.phase_z).sin()
924            }
925            2 => {
926                -self.freq_z
927                    * self.freq_z
928                    * (self.freq_x * x + self.phase_x).sin()
929                    * (self.freq_y * y + self.phase_y).cos()
930                    * (self.freq_z * z + self.phase_z).sin()
931            }
932            _ => F::zero(),
933        }
934    }
935
936    fn dimension(&self) -> usize {
937        3
938    }
939}
940
941/// System of equations verification (simplified without trait objects)
942#[derive(Debug, Clone)]
943pub struct SystemVerification<F: IntegrateFloat> {
944    /// Number of components in the system
945    pub system_size: usize,
946    /// Component names for identification
947    pub component_names: Vec<String>,
948    /// Phantom data to maintain type parameter
949    _phantom: std::marker::PhantomData<F>,
950}
951
952impl<F: IntegrateFloat> SystemVerification<F> {
953    /// Create new system verification
954    pub fn new(systemsize: usize) -> Self {
955        let component_names = (0..systemsize).map(|i| format!("Component {i}")).collect();
956        Self {
957            system_size: systemsize,
958            component_names,
959            _phantom: std::marker::PhantomData,
960        }
961    }
962
963    /// Create with custom component names
964    pub fn with_names(componentnames: Vec<String>) -> Self {
965        let system_size = componentnames.len();
966        Self {
967            system_size,
968            component_names: componentnames,
969            _phantom: std::marker::PhantomData,
970        }
971    }
972
973    /// Verify system solution using provided exact solutions
974    pub fn verify_system<S1, S2>(
975        &self,
976        exact_solutions: &[S1],
977        numerical_solutions: &[S2],
978        coordinates: &[F],
979    ) -> Vec<F>
980    where
981        S1: ExactSolution<F>,
982        S2: Fn(&[F]) -> F,
983    {
984        assert_eq!(exact_solutions.len(), self.system_size);
985        assert_eq!(numerical_solutions.len(), self.system_size);
986
987        let mut errors = Vec::with_capacity(self.system_size);
988        for i in 0..self.system_size {
989            let exact = exact_solutions[i].evaluate(coordinates);
990            let numerical = numerical_solutions[i](coordinates);
991            errors.push((exact - numerical).abs());
992        }
993        errors
994    }
995}
996
997/// Automated verification workflow for numerical methods
998pub struct VerificationWorkflow<F: IntegrateFloat> {
999    /// Test cases to verify
1000    pub test_cases: Vec<VerificationTestCase<F>>,
1001}
1002
1003#[derive(Debug, Clone)]
1004pub struct VerificationTestCase<F: IntegrateFloat> {
1005    /// Test name
1006    pub name: String,
1007    /// Expected order of accuracy
1008    pub expected_order: F,
1009    /// Tolerance for order verification
1010    pub order_tolerance: F,
1011    /// Grid sizes to test
1012    pub grid_sizes: Vec<F>,
1013    /// Expected errors (if known)
1014    pub expected_errors: Option<Vec<F>>,
1015}
1016
1017impl<F: IntegrateFloat> VerificationWorkflow<F> {
1018    /// Create new verification workflow
1019    pub fn new() -> Self {
1020        Self {
1021            test_cases: Vec::new(),
1022        }
1023    }
1024
1025    /// Add test case to workflow
1026    pub fn add_test_case(&mut self, testcase: VerificationTestCase<F>) {
1027        self.test_cases.push(testcase);
1028    }
1029
1030    /// Run all verification tests
1031    pub fn run_verification<S: Fn(&[F]) -> IntegrateResult<F>>(
1032        &self,
1033        solver: S,
1034    ) -> Vec<VerificationResult<F>> {
1035        let mut results = Vec::new();
1036
1037        for test_case in &self.test_cases {
1038            let mut errors = Vec::new();
1039
1040            for &grid_size in &test_case.grid_sizes {
1041                match solver(&[grid_size]) {
1042                    Ok(error) => errors.push(error),
1043                    Err(_) => {
1044                        results.push(VerificationResult {
1045                            test_name: test_case.name.clone(),
1046                            passed: false,
1047                            computed_order: None,
1048                            error_message: Some("Solver failed".to_string()),
1049                        });
1050                        break;
1051                    }
1052                }
1053            }
1054
1055            if errors.len() == test_case.grid_sizes.len() {
1056                match ConvergenceAnalysis::compute_order(test_case.grid_sizes.clone(), errors) {
1057                    Ok(analysis) => {
1058                        let passed = analysis
1059                            .verify_order(test_case.expected_order, test_case.order_tolerance);
1060                        results.push(VerificationResult {
1061                            test_name: test_case.name.clone(),
1062                            passed,
1063                            computed_order: Some(analysis.order),
1064                            error_message: if passed {
1065                                None
1066                            } else {
1067                                Some("Order verification failed".to_string())
1068                            },
1069                        });
1070                    }
1071                    Err(e) => {
1072                        results.push(VerificationResult {
1073                            test_name: test_case.name.clone(),
1074                            passed: false,
1075                            computed_order: None,
1076                            error_message: Some(format!("Convergence analysis failed: {e}")),
1077                        });
1078                    }
1079                }
1080            }
1081        }
1082
1083        results
1084    }
1085}
1086
1087#[derive(Debug, Clone)]
1088pub struct VerificationResult<F: IntegrateFloat> {
1089    /// Test name
1090    pub test_name: String,
1091    /// Whether test passed
1092    pub passed: bool,
1093    /// Computed order of accuracy
1094    pub computed_order: Option<F>,
1095    /// Error message if test failed
1096    pub error_message: Option<String>,
1097}
1098
1099impl<F: IntegrateFloat> Default for VerificationWorkflow<F> {
1100    fn default() -> Self {
1101        Self::new()
1102    }
1103}
1104
1105#[allow(dead_code)]
1106pub fn exponential_solution<F: IntegrateFloat>(
1107    amplitude: F,
1108    decay_rate: F,
1109) -> ExponentialSolution<F> {
1110    ExponentialSolution::simple(amplitude, decay_rate)
1111}
1112
1113#[allow(dead_code)]
1114pub fn trigonometric_solution_3d<F: IntegrateFloat>(
1115    freq_x: F,
1116    freq_y: F,
1117    freq_z: F,
1118) -> TrigonometricSolution3D<F> {
1119    TrigonometricSolution3D::simple(freq_x, freq_y, freq_z)
1120}
1121
1122#[allow(dead_code)]
1123pub fn combined_solution<F: IntegrateFloat>(dimension: usize) -> CombinedSolution<F> {
1124    CombinedSolution::new(dimension)
1125}
1126
1127#[cfg(test)]
1128mod tests {
1129    use super::*;
1130    use approx::assert_abs_diff_eq;
1131    use scirs2_core::ndarray::Array1;
1132
1133    #[test]
1134    fn test_polynomial_solution() {
1135        // Test polynomial: 1 + 2t + 3t²
1136        let poly = polynomial_solution(vec![1.0, 2.0, 3.0]);
1137
1138        // Test evaluation
1139        assert_abs_diff_eq!(poly.evaluate(&[0.0]), 1.0);
1140        assert_abs_diff_eq!(poly.evaluate(&[1.0]), 6.0); // 1 + 2 + 3
1141
1142        // Test derivative: 2 + 6t
1143        assert_abs_diff_eq!(poly.derivative(&[0.0], 0), 2.0);
1144        assert_abs_diff_eq!(poly.derivative(&[1.0], 0), 8.0);
1145
1146        // Test second derivative: 6
1147        assert_abs_diff_eq!(poly.second_derivative(&[0.0], 0), 6.0);
1148        assert_abs_diff_eq!(poly.second_derivative(&[1.0], 0), 6.0);
1149    }
1150
1151    #[test]
1152    fn test_trigonometric_solution_2d() {
1153        // Test sin(x) * cos(y)
1154        let trig = trigonometric_solution_2d(1.0, 1.0);
1155
1156        // Test evaluation
1157        assert_abs_diff_eq!(trig.evaluate(&[0.0, 0.0]), 0.0); // sin(0) * cos(0) = 0 * 1 = 0
1158        assert_abs_diff_eq!(trig.evaluate(&[PI / 2.0, 0.0]), 1.0); // sin(π/2) * cos(0) = 1 * 1 = 1
1159
1160        // Test derivatives
1161        // ∂/∂x[sin(x)cos(y)] = cos(x)cos(y)
1162        assert_abs_diff_eq!(trig.derivative(&[0.0, 0.0], 0), 1.0); // cos(0)*cos(0) = 1
1163
1164        // ∂/∂y[sin(x)cos(y)] = -sin(x)sin(y)
1165        assert_abs_diff_eq!(trig.derivative(&[0.0, 0.0], 1), 0.0); // -sin(0)*sin(0) = 0
1166    }
1167
1168    #[test]
1169    fn test_mms_ode_problem() {
1170        // Test with y = t²
1171        let poly = polynomial_solution(vec![0.0, 0.0, 1.0]); // t²
1172        let problem = MMSODEProblem::new(poly, [0.0, 1.0]);
1173
1174        // Initial condition should be 0
1175        assert_abs_diff_eq!(problem.initial_condition(), 0.0);
1176
1177        // Source term should be y' = 2t
1178        assert_abs_diff_eq!(problem.source_term(0.0), 0.0);
1179        assert_abs_diff_eq!(problem.source_term(1.0), 2.0);
1180
1181        // Exact solution at t=0.5 should be 0.25
1182        assert_abs_diff_eq!(problem.exact_at(0.5), 0.25);
1183    }
1184
1185    #[test]
1186    fn test_convergence_analysis() {
1187        // Test with theoretical O(h²) convergence
1188        let grid_sizes = vec![0.1, 0.05, 0.025, 0.0125];
1189        let errors = vec![0.01, 0.0025, 0.000625, 0.00015625]; // ∝ h²
1190
1191        let analysis = ConvergenceAnalysis::compute_order(grid_sizes, errors).unwrap();
1192
1193        // Should detect order ≈ 2
1194        assert!((analysis.order - 2.0_f64).abs() < 0.1);
1195        assert!(analysis.verify_order(2.0, 0.2));
1196    }
1197
1198    #[test]
1199    fn test_error_analysis() {
1200        let exact = Array1::from_vec(vec![1.0, 2.0, 3.0, 4.0]);
1201        let numerical = Array1::from_vec(vec![1.1, 1.9, 3.05, 3.98]);
1202
1203        let l2_error = ErrorAnalysis::l2_norm(exact.view(), numerical.view()).unwrap();
1204        let max_error = ErrorAnalysis::max_norm(exact.view(), numerical.view()).unwrap();
1205
1206        // L2 error should be around sqrt((0.1² + 0.1² + 0.05² + 0.02²)/4)
1207        assert!(l2_error > 0.0 && l2_error < 0.2);
1208
1209        // Max error should be 0.1
1210        assert_abs_diff_eq!(max_error, 0.1, epsilon = 1e-10);
1211    }
1212
1213    #[test]
1214    fn test_exponential_solution() {
1215        // Test exponential solution: 2 * exp(-3t)
1216        let exp_sol = exponential_solution(2.0, -3.0);
1217
1218        // Test evaluation
1219        assert_abs_diff_eq!(exp_sol.evaluate(&[0.0]), 2.0);
1220        assert_abs_diff_eq!(exp_sol.evaluate(&[1.0]), 2.0 * (-3.0_f64).exp());
1221
1222        // Test derivative: d/dt[2 * exp(-3t)] = -6 * exp(-3t)
1223        assert_abs_diff_eq!(exp_sol.derivative(&[0.0], 0), -6.0);
1224        assert_abs_diff_eq!(exp_sol.derivative(&[1.0], 0), -6.0 * (-3.0_f64).exp());
1225
1226        // Test second derivative: d²/dt²[2 * exp(-3t)] = 18 * exp(-3t)
1227        assert_abs_diff_eq!(exp_sol.second_derivative(&[0.0], 0), 18.0);
1228        assert_abs_diff_eq!(
1229            exp_sol.second_derivative(&[1.0], 0),
1230            18.0 * (-3.0_f64).exp()
1231        );
1232    }
1233
1234    #[test]
1235    fn test_combined_solution() {
1236        // Test combined solution: polynomial + exponential
1237        let poly = polynomial_solution(vec![1.0, 2.0]); // 1 + 2t
1238        let exp = exponential_solution(1.0, -1.0); // exp(-t)
1239
1240        let combined = combined_solution(1)
1241            .with_polynomial(poly)
1242            .with_exponential(exp);
1243
1244        // At t=0: (1 + 2*0) + exp(-0) = 1 + 1 = 2
1245        assert_abs_diff_eq!(combined.evaluate(&[0.0]), 2.0);
1246
1247        // At t=1: (1 + 2*1) + exp(-1) = 3 + exp(-1)
1248        let expected = 3.0 + (-1.0_f64).exp();
1249        assert_abs_diff_eq!(combined.evaluate(&[1.0]), expected, epsilon = 1e-10);
1250
1251        // Test derivative at t=0: d/dt[1 + 2t + exp(-t)] = 2 - exp(-t)
1252        // At t=0: 2 - 1 = 1
1253        assert_abs_diff_eq!(combined.derivative(&[0.0], 0), 1.0);
1254    }
1255
1256    #[test]
1257    fn test_trigonometric_solution_3d() {
1258        // Test sin(x) * cos(y) * sin(z)
1259        let trig3d = trigonometric_solution_3d(1.0, 1.0, 1.0);
1260
1261        // Test evaluation at (π/2, 0, π/2)
1262        // sin(π/2) * cos(0) * sin(π/2) = 1 * 1 * 1 = 1
1263        assert_abs_diff_eq!(
1264            trig3d.evaluate(&[PI / 2.0, 0.0, PI / 2.0]),
1265            1.0,
1266            epsilon = 1e-10
1267        );
1268
1269        // Test derivative with respect to x at (0, 0, π/2)
1270        // ∂/∂x[sin(x)cos(y)sin(z)] = cos(x)cos(y)sin(z)
1271        // cos(0)*cos(0)*sin(π/2) = 1*1*1 = 1
1272        assert_abs_diff_eq!(
1273            trig3d.derivative(&[0.0, 0.0, PI / 2.0], 0),
1274            1.0,
1275            epsilon = 1e-10
1276        );
1277
1278        // Test second derivative with respect to x at (0, 0, π/2)
1279        // ∂²/∂x²[sin(x)cos(y)sin(z)] = -sin(x)cos(y)sin(z)
1280        // -sin(0)*cos(0)*sin(π/2) = 0
1281        assert_abs_diff_eq!(
1282            trig3d.second_derivative(&[0.0, 0.0, PI / 2.0], 0),
1283            0.0,
1284            epsilon = 1e-10
1285        );
1286    }
1287
1288    #[test]
1289    fn test_3d_poisson_problem() {
1290        // Test 3D Poisson problem with trigonometric solution
1291        let exact = trigonometric_solution_3d(PI, PI, PI);
1292        let problem = MMSPDEProblem::new_poisson_3d(exact, [0.0, 1.0], [0.0, 1.0], [0.0, 1.0]);
1293
1294        // Test that it's recognized as 3D
1295        assert!(problem.is_3d());
1296
1297        // Test domain
1298        let (dx, dy, dz) = problem.domain_3d();
1299        assert_eq!(dx, [0.0, 1.0]);
1300        assert_eq!(dy, [0.0, 1.0]);
1301        assert_eq!(dz, [0.0, 1.0]);
1302
1303        // Test source term computation
1304        // For -∇²u = f where u = sin(πx)cos(πy)sin(πz)
1305        // ∇²u = -π²sin(πx)cos(πy)sin(πz) - π²sin(πx)cos(πy)sin(πz) - π²sin(πx)cos(πy)sin(πz)
1306        //     = -3π²sin(πx)cos(πy)sin(πz) = -3π²u
1307        // So f = -∇²u = 3π²u
1308        let coords = [0.5, 0.0, 0.5]; // Where cos(πy) = cos(0) = 1
1309        let u_exact = problem.exact_at_3d(coords[0], coords[1], coords[2]);
1310        let f_computed = problem.source_term(&coords);
1311        let f_expected = 3.0 * PI * PI * u_exact;
1312        assert_abs_diff_eq!(f_computed, f_expected, epsilon = 1e-10);
1313    }
1314
1315    #[test]
1316    fn test_helmholtz_problem() {
1317        // Test 2D Helmholtz problem: ∇²u + k²u = f
1318        let exact = trigonometric_solution_2d(PI, PI);
1319        let k = 2.0;
1320        let problem = MMSPDEProblem::new_helmholtz_2d(exact, [0.0, 1.0], [0.0, 1.0], k);
1321
1322        // For u = sin(πx)cos(πy), ∇²u = -2π²sin(πx)cos(πy) = -2π²u
1323        // So f = ∇²u + k²u = -2π²u + k²u = (k² - 2π²)u
1324        let coords = [0.5, 0.0]; // Where cos(πy) = cos(0) = 1
1325        let u_exact = problem.exact_at(coords[0], coords[1]);
1326        let f_computed = problem.source_term(&coords);
1327        let f_expected = (k * k - 2.0 * PI * PI) * u_exact;
1328        assert_abs_diff_eq!(f_computed, f_expected, epsilon = 1e-10);
1329    }
1330
1331    #[test]
1332    fn test_verification_workflow() {
1333        let mut workflow = VerificationWorkflow::new();
1334
1335        // Add a test case for second-order method
1336        let test_case = VerificationTestCase {
1337            name: "Second-order test".to_string(),
1338            expected_order: 2.0,
1339            order_tolerance: 0.1,
1340            grid_sizes: vec![0.1, 0.05, 0.025],
1341            expected_errors: None,
1342        };
1343        workflow.add_test_case(test_case);
1344
1345        // Mock solver that produces second-order errors
1346        let mock_solver = |grid_sizes: &[f64]| -> IntegrateResult<f64> {
1347            let h = grid_sizes[0];
1348            Ok(0.1 * h * h) // O(h²) error
1349        };
1350
1351        let results = workflow.run_verification(mock_solver);
1352        assert_eq!(results.len(), 1);
1353        assert!(results[0].passed);
1354        assert!(results[0].computed_order.unwrap() > 1.8);
1355        assert!(results[0].computed_order.unwrap() < 2.2);
1356    }
1357
1358    #[test]
1359    fn test_system_verification() {
1360        // Test system with two components: polynomial and trigonometric
1361        let poly = polynomial_solution(vec![1.0, 2.0]);
1362        let trig = trigonometric_solution_2d(1.0, 1.0);
1363
1364        let system = SystemVerification::<f64>::new(2);
1365        assert_eq!(system.system_size, 2);
1366
1367        // Test verification at (0, 0) - handle solutions separately due to different types
1368        let coords = vec![0.0, 0.0];
1369
1370        // Test polynomial solution
1371        let poly_exact = poly.evaluate(&coords);
1372        let poly_numerical = 1.0 + 2.0 * coords[0]; // approximate polynomial
1373
1374        // Test trigonometric solution
1375        let trig_exact = trig.evaluate(&coords);
1376        let trig_numerical = coords[0] * coords[1]; // approximate trigonometric
1377
1378        // Calculate errors manually since we can't use mixed-type arrays
1379        let poly_error = (poly_exact as f64 - poly_numerical).abs() as f64;
1380        let trig_error = (trig_exact as f64 - trig_numerical).abs() as f64;
1381
1382        // At (0,0): exact poly = 1, numerical = 1, error = 0
1383        assert_abs_diff_eq!(poly_error, 0.0);
1384        // At (0,0): exact trig = 0, numerical = 0, error = 0
1385        assert_abs_diff_eq!(trig_error, 0.0);
1386
1387        // Test with custom names
1388        let named_system: SystemVerification<f64> = SystemVerification::with_names(vec![
1389            "Polynomial".to_string(),
1390            "Trigonometric".to_string(),
1391        ]);
1392        assert_eq!(named_system.component_names[0], "Polynomial");
1393        assert_eq!(named_system.component_names[1], "Trigonometric");
1394    }
1395}
1396
1397/// Advanced verification framework with automatic error estimation
1398pub struct AdvancedVerificationFramework {
1399    /// Grid refinement strategy
1400    pub refinement_strategy: RefinementStrategy,
1401    /// Error estimation method
1402    pub error_estimation_method: ErrorEstimationMethod,
1403    /// Convergence criteria
1404    pub convergence_criteria: ConvergenceCriteria,
1405    /// Statistical analysis settings
1406    pub statistical_analysis: bool,
1407}
1408
1409/// Grid refinement strategies
1410#[derive(Debug, Clone, Copy)]
1411pub enum RefinementStrategy {
1412    /// Uniform refinement (h → h/2)
1413    Uniform,
1414    /// Adaptive refinement based on error indicators
1415    Adaptive,
1416    /// Custom refinement ratios
1417    Custom(f64),
1418}
1419
1420/// Error estimation methods
1421#[derive(Debug, Clone, Copy)]
1422pub enum ErrorEstimationMethod {
1423    /// Richardson extrapolation
1424    Richardson,
1425    /// Embedded method comparison
1426    Embedded,
1427    /// Bootstrap sampling
1428    Bootstrap,
1429    /// Cross-validation
1430    CrossValidation,
1431}
1432
1433/// Convergence criteria for verification
1434#[derive(Debug, Clone)]
1435pub struct ConvergenceCriteria {
1436    /// Minimum number of grid levels
1437    pub min_levels: usize,
1438    /// Maximum number of grid levels
1439    pub max_levels: usize,
1440    /// Required R² for order estimation
1441    pub min_r_squared: f64,
1442    /// Tolerance for order estimation
1443    pub order_tolerance: f64,
1444    /// Target accuracy
1445    pub target_accuracy: f64,
1446}
1447
1448impl Default for ConvergenceCriteria {
1449    fn default() -> Self {
1450        Self {
1451            min_levels: 3,
1452            max_levels: 8,
1453            min_r_squared: 0.95,
1454            order_tolerance: 0.2,
1455            target_accuracy: 1e-10,
1456        }
1457    }
1458}
1459
1460impl Default for AdvancedVerificationFramework {
1461    fn default() -> Self {
1462        Self {
1463            refinement_strategy: RefinementStrategy::Uniform,
1464            error_estimation_method: ErrorEstimationMethod::Richardson,
1465            convergence_criteria: ConvergenceCriteria::default(),
1466            statistical_analysis: true,
1467        }
1468    }
1469}