Skip to main content

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