bem/core/
bem_solver.rs

1//! High-level BEM Solver API
2//!
3//! This module provides a unified, high-level interface for acoustic BEM simulations.
4//! It integrates mesh generation, system assembly, linear solving, and post-processing.
5//!
6//! The Direct solver method works in both native and WASM modes. When the `native` feature
7//! is enabled, it uses optimized BLAS/LAPACK. Otherwise, it falls back to a pure Rust
8//! LU factorization implementation.
9//!
10//! # Example
11//!
12//! ```ignore
13//! use bem::core::{BemSolver, BemProblem, IncidentField};
14//!
15//! // Create a rigid sphere scattering problem
16//! let problem = BemProblem::rigid_sphere_scattering(
17//!     0.1,        // radius
18//!     1000.0,     // frequency
19//!     343.0,      // speed of sound
20//!     1.21,       // density
21//! );
22//!
23//! // Configure solver
24//! let solver = BemSolver::new()
25//!     .with_mesh_refinement(3)
26//!     .with_solver_method(SolverMethod::Direct);
27//!
28//! // Solve
29//! let solution = solver.solve(&problem)?;
30//!
31//! // Evaluate field at a point
32//! let pressure = solution.evaluate_pressure(&[0.0, 0.0, 0.2]);
33//! ```
34
35use ndarray::{Array1, Array2};
36use num_complex::Complex64;
37use std::f64::consts::PI;
38
39use crate::core::assembly::tbem::build_tbem_system_with_beta;
40use crate::core::incident::IncidentField;
41use crate::core::mesh::generators::{generate_icosphere_mesh, generate_sphere_mesh};
42use crate::core::postprocess::pressure::{FieldPoint, compute_total_field};
43use crate::core::types::{BoundaryCondition, Element, Mesh, PhysicsParams};
44use solvers::direct::lu_solve;
45
46/// Solver method for the linear system
47#[derive(Debug, Clone, Copy, PartialEq, Eq, Default)]
48pub enum SolverMethod {
49    /// Direct LU factorization (best for small problems)
50    #[default]
51    Direct,
52    /// Conjugate Gradient Squared (iterative)
53    Cgs,
54    /// BiCGSTAB (iterative, more stable than CGS)
55    BiCgStab,
56}
57
58/// Assembly method for the BEM matrix
59#[derive(Debug, Clone, Copy, PartialEq, Eq, Default)]
60pub enum AssemblyMethod {
61    /// Traditional BEM with O(N^2) dense matrix
62    #[default]
63    Tbem,
64    /// Single-Level Fast Multipole Method
65    Slfmm,
66    /// Multi-Level Fast Multipole Method
67    Mlfmm,
68}
69
70/// Boundary condition type for the problem
71#[derive(Debug, Clone, Copy, PartialEq, Eq)]
72pub enum BoundaryConditionType {
73    /// Rigid surface (zero normal velocity)
74    Rigid,
75    /// Soft surface (zero pressure)
76    Soft,
77    /// Impedance boundary condition
78    Impedance,
79}
80
81/// Definition of a BEM problem
82#[derive(Debug, Clone)]
83pub struct BemProblem {
84    /// Problem geometry mesh
85    pub mesh: Mesh,
86    /// Physical parameters
87    pub physics: PhysicsParams,
88    /// Incident field
89    pub incident_field: IncidentField,
90    /// Boundary condition type
91    pub bc_type: BoundaryConditionType,
92    /// Use Burton-Miller formulation (recommended for exterior problems)
93    pub use_burton_miller: bool,
94}
95
96impl BemProblem {
97    /// Create a rigid sphere scattering problem with plane wave incidence
98    ///
99    /// # Arguments
100    /// * `radius` - Sphere radius (m)
101    /// * `frequency` - Excitation frequency (Hz)
102    /// * `speed_of_sound` - Speed of sound (m/s)
103    /// * `density` - Medium density (kg/m�)
104    pub fn rigid_sphere_scattering(
105        radius: f64,
106        frequency: f64,
107        speed_of_sound: f64,
108        density: f64,
109    ) -> Self {
110        // Determine mesh resolution based on ka
111        let k = 2.0 * PI * frequency / speed_of_sound;
112        let ka = k * radius;
113
114        // Rule of thumb: ~10 elements per wavelength
115        // For a sphere, this translates to subdivision level
116        let subdivisions = if ka < 1.0 {
117            2 // Low frequency: coarse mesh ok
118        } else if ka < 5.0 {
119            3 // Medium frequency
120        } else {
121            4 // High frequency: need finer mesh
122        };
123
124        let mesh = generate_icosphere_mesh(radius, subdivisions);
125        let physics = PhysicsParams::new(frequency, speed_of_sound, density, false);
126        let incident_field = IncidentField::plane_wave_z();
127
128        Self {
129            mesh,
130            physics,
131            incident_field,
132            bc_type: BoundaryConditionType::Rigid,
133            use_burton_miller: true,
134        }
135    }
136
137    /// Create a rigid sphere scattering problem with custom mesh resolution
138    pub fn rigid_sphere_scattering_custom(
139        radius: f64,
140        frequency: f64,
141        speed_of_sound: f64,
142        density: f64,
143        n_theta: usize,
144        n_phi: usize,
145    ) -> Self {
146        let mesh = generate_sphere_mesh(radius, n_theta, n_phi);
147        let physics = PhysicsParams::new(frequency, speed_of_sound, density, false);
148        let incident_field = IncidentField::plane_wave_z();
149
150        Self {
151            mesh,
152            physics,
153            incident_field,
154            bc_type: BoundaryConditionType::Rigid,
155            use_burton_miller: true,
156        }
157    }
158
159    /// Set the incident field
160    pub fn with_incident_field(mut self, field: IncidentField) -> Self {
161        self.incident_field = field;
162        self
163    }
164
165    /// Set the boundary condition type
166    pub fn with_boundary_condition(mut self, bc_type: BoundaryConditionType) -> Self {
167        self.bc_type = bc_type;
168        self
169    }
170
171    /// Enable/disable Burton-Miller formulation
172    pub fn with_burton_miller(mut self, use_bm: bool) -> Self {
173        self.use_burton_miller = use_bm;
174        self
175    }
176
177    /// Get the wave number times radius (ka)
178    pub fn ka(&self) -> f64 {
179        self.physics.wave_number * self.mesh_radius()
180    }
181
182    /// Estimate mesh radius from bounding box
183    fn mesh_radius(&self) -> f64 {
184        // Estimate radius from mesh nodes
185        let mut max_r = 0.0f64;
186        for i in 0..self.mesh.nodes.nrows() {
187            let r = (self.mesh.nodes[[i, 0]].powi(2)
188                + self.mesh.nodes[[i, 1]].powi(2)
189                + self.mesh.nodes[[i, 2]].powi(2))
190            .sqrt();
191            max_r = max_r.max(r);
192        }
193        max_r
194    }
195}
196
197/// BEM solver configuration
198#[derive(Debug, Clone)]
199pub struct BemSolver {
200    /// Linear solver method
201    pub solver_method: SolverMethod,
202    /// Matrix assembly method
203    pub assembly_method: AssemblyMethod,
204    /// Maximum iterations for iterative solvers
205    pub max_iterations: usize,
206    /// Tolerance for iterative solvers
207    pub tolerance: f64,
208    /// Verbose output
209    pub verbose: bool,
210    /// Burton-Miller β scale factor (default: 4.0 for best accuracy at ka ~ 1)
211    pub beta_scale: f64,
212}
213
214impl Default for BemSolver {
215    fn default() -> Self {
216        Self {
217            solver_method: SolverMethod::Direct,
218            assembly_method: AssemblyMethod::Tbem,
219            max_iterations: 1000,
220            tolerance: 1e-8,
221            verbose: false,
222            beta_scale: 4.0, // Empirically optimal for ka ~ 1
223        }
224    }
225}
226
227impl BemSolver {
228    /// Create a new solver with default settings
229    pub fn new() -> Self {
230        Self::default()
231    }
232
233    /// Set the linear solver method
234    pub fn with_solver_method(mut self, method: SolverMethod) -> Self {
235        self.solver_method = method;
236        self
237    }
238
239    /// Set the assembly method
240    pub fn with_assembly_method(mut self, method: AssemblyMethod) -> Self {
241        self.assembly_method = method;
242        self
243    }
244
245    /// Set maximum iterations for iterative solvers
246    pub fn with_max_iterations(mut self, max_iter: usize) -> Self {
247        self.max_iterations = max_iter;
248        self
249    }
250
251    /// Set tolerance for iterative solvers
252    pub fn with_tolerance(mut self, tol: f64) -> Self {
253        self.tolerance = tol;
254        self
255    }
256
257    /// Enable verbose output
258    pub fn with_verbose(mut self, verbose: bool) -> Self {
259        self.verbose = verbose;
260        self
261    }
262
263    /// Solve a BEM problem
264    ///
265    /// # Arguments
266    /// * `problem` - The BEM problem definition
267    ///
268    /// # Returns
269    /// Solution containing surface pressures and methods to evaluate fields
270    pub fn solve(&self, problem: &BemProblem) -> Result<BemSolution, BemError> {
271        if self.verbose {
272            log::info!(
273                "Solving BEM problem: {} elements, ka = {:.3}",
274                problem.mesh.elements.len(),
275                problem.ka()
276            );
277        }
278
279        // Step 1: Prepare elements with boundary conditions
280        let elements = self.prepare_elements(problem);
281
282        // Step 2: Assemble system matrix and RHS
283        let (matrix, rhs) =
284            self.assemble_system(&elements, &problem.mesh.nodes, &problem.physics)?;
285
286        // Step 3: Add incident field contribution to RHS
287        let rhs = self.add_incident_field_rhs(
288            rhs,
289            &elements,
290            &problem.incident_field,
291            &problem.physics,
292            problem.use_burton_miller,
293        );
294
295        // Step 4: Solve linear system
296        let surface_pressure = self.solve_linear_system(&matrix, &rhs)?;
297
298        if self.verbose {
299            log::info!(
300                "Solution complete. Max surface pressure: {:.6}",
301                surface_pressure
302                    .iter()
303                    .map(|p| p.norm())
304                    .fold(0.0f64, f64::max)
305            );
306        }
307
308        Ok(BemSolution {
309            surface_pressure,
310            elements,
311            nodes: problem.mesh.nodes.clone(),
312            incident_field: problem.incident_field.clone(),
313            physics: problem.physics.clone(),
314        })
315    }
316
317    /// Prepare elements with boundary conditions
318    fn prepare_elements(&self, problem: &BemProblem) -> Vec<Element> {
319        let mut elements = problem.mesh.elements.clone();
320
321        // Set boundary conditions based on problem type
322        let bc = match problem.bc_type {
323            BoundaryConditionType::Rigid => {
324                // Zero normal velocity (Neumann BC)
325                BoundaryCondition::Velocity(vec![Complex64::new(0.0, 0.0)])
326            }
327            BoundaryConditionType::Soft => {
328                // Zero pressure (Dirichlet BC)
329                BoundaryCondition::Pressure(vec![Complex64::new(0.0, 0.0)])
330            }
331            BoundaryConditionType::Impedance => {
332                // Default impedance (plane wave)
333                let z0 = problem.physics.density * problem.physics.speed_of_sound;
334                BoundaryCondition::VelocityWithAdmittance {
335                    velocity: vec![Complex64::new(0.0, 0.0)],
336                    admittance: Complex64::new(1.0 / z0, 0.0),
337                }
338            }
339        };
340
341        // Assign BC and DOF addresses to each element
342        for (i, elem) in elements.iter_mut().enumerate() {
343            elem.boundary_condition = bc.clone();
344            elem.dof_addresses = vec![i];
345        }
346
347        elements
348    }
349
350    /// Assemble the BEM system matrix
351    fn assemble_system(
352        &self,
353        elements: &[Element],
354        nodes: &Array2<f64>,
355        physics: &PhysicsParams,
356    ) -> Result<(Array2<Complex64>, Array1<Complex64>), BemError> {
357        match self.assembly_method {
358            AssemblyMethod::Tbem => {
359                // Use scaled Burton-Miller β for better accuracy
360                let beta = physics.burton_miller_beta_scaled(self.beta_scale);
361                let system = build_tbem_system_with_beta(elements, nodes, physics, beta);
362                Ok((system.matrix, system.rhs))
363            }
364            AssemblyMethod::Slfmm | AssemblyMethod::Mlfmm => {
365                // FMM methods not yet integrated into high-level API
366                Err(BemError::NotImplemented(
367                    "FMM assembly not yet available in high-level API".to_string(),
368                ))
369            }
370        }
371    }
372
373    /// Add incident field contribution to RHS
374    fn add_incident_field_rhs(
375        &self,
376        mut rhs: Array1<Complex64>,
377        elements: &[Element],
378        incident_field: &IncidentField,
379        physics: &PhysicsParams,
380        use_burton_miller: bool,
381    ) -> Array1<Complex64> {
382        // Collect element centers and normals
383        let n = elements.len();
384        let mut centers = Array2::zeros((n, 3));
385        let mut normals = Array2::zeros((n, 3));
386
387        for (i, elem) in elements.iter().enumerate() {
388            for j in 0..3 {
389                centers[[i, j]] = elem.center[j];
390                normals[[i, j]] = elem.normal[j];
391            }
392        }
393
394        // Compute incident field RHS contribution with scaled β
395        let incident_rhs = if use_burton_miller {
396            let beta = physics.burton_miller_beta_scaled(self.beta_scale);
397            incident_field.compute_rhs_with_beta(&centers, &normals, physics, beta)
398        } else {
399            incident_field.compute_rhs(&centers, &normals, physics, false)
400        };
401
402        // Add to system RHS
403        rhs = rhs + incident_rhs;
404
405        rhs
406    }
407
408    /// Solve the linear system
409    fn solve_linear_system(
410        &self,
411        matrix: &Array2<Complex64>,
412        rhs: &Array1<Complex64>,
413    ) -> Result<Array1<Complex64>, BemError> {
414        match self.solver_method {
415            SolverMethod::Direct => {
416                // Uses BLAS when native feature is enabled, pure Rust fallback otherwise
417                lu_solve(matrix, rhs).map_err(|e| BemError::SolverFailed(e.to_string()))
418            }
419            SolverMethod::Cgs | SolverMethod::BiCgStab => {
420                // Iterative solvers not yet integrated
421                Err(BemError::NotImplemented(
422                    "Iterative solvers not yet available in high-level API".to_string(),
423                ))
424            }
425        }
426    }
427}
428
429/// Solution of a BEM problem
430#[derive(Debug, Clone)]
431pub struct BemSolution {
432    /// Surface pressure at each element
433    pub surface_pressure: Array1<Complex64>,
434    /// Elements used in the solution
435    pub elements: Vec<Element>,
436    /// Node coordinates
437    pub nodes: Array2<f64>,
438    /// Incident field used
439    pub incident_field: IncidentField,
440    /// Physics parameters
441    pub physics: PhysicsParams,
442}
443
444impl BemSolution {
445    /// Evaluate total pressure at a single point
446    pub fn evaluate_pressure(&self, point: &[f64; 3]) -> Complex64 {
447        let eval_points =
448            Array2::from_shape_vec((1, 3), vec![point[0], point[1], point[2]]).unwrap();
449
450        let field_points = compute_total_field(
451            &eval_points,
452            &self.elements,
453            &self.nodes,
454            &self.surface_pressure,
455            None,
456            &self.incident_field,
457            &self.physics,
458        );
459
460        field_points[0].p_total
461    }
462
463    /// Evaluate total pressure at multiple points
464    pub fn evaluate_pressure_field(&self, points: &Array2<f64>) -> Vec<FieldPoint> {
465        compute_total_field(
466            points,
467            &self.elements,
468            &self.nodes,
469            &self.surface_pressure,
470            None,
471            &self.incident_field,
472            &self.physics,
473        )
474    }
475
476    /// Get max surface pressure magnitude
477    pub fn max_surface_pressure(&self) -> f64 {
478        self.surface_pressure
479            .iter()
480            .map(|p| p.norm())
481            .fold(0.0f64, f64::max)
482    }
483
484    /// Get mean surface pressure magnitude
485    pub fn mean_surface_pressure(&self) -> f64 {
486        let sum: f64 = self.surface_pressure.iter().map(|p| p.norm()).sum();
487        sum / self.surface_pressure.len() as f64
488    }
489
490    /// Number of DOFs in the solution
491    pub fn num_dofs(&self) -> usize {
492        self.surface_pressure.len()
493    }
494}
495
496/// BEM solver errors
497#[derive(Debug, Clone)]
498pub enum BemError {
499    /// Feature not yet implemented
500    NotImplemented(String),
501    /// Linear solver failed
502    SolverFailed(String),
503    /// Invalid mesh
504    InvalidMesh(String),
505    /// Invalid parameters
506    InvalidParameters(String),
507}
508
509impl std::fmt::Display for BemError {
510    fn fmt(&self, f: &mut std::fmt::Formatter<'_>) -> std::fmt::Result {
511        match self {
512            BemError::NotImplemented(msg) => write!(f, "Not implemented: {}", msg),
513            BemError::SolverFailed(msg) => write!(f, "Solver failed: {}", msg),
514            BemError::InvalidMesh(msg) => write!(f, "Invalid mesh: {}", msg),
515            BemError::InvalidParameters(msg) => write!(f, "Invalid parameters: {}", msg),
516        }
517    }
518}
519
520impl std::error::Error for BemError {}
521
522#[cfg(test)]
523mod tests {
524    use super::*;
525
526    #[test]
527    fn test_bem_problem_creation() {
528        let problem = BemProblem::rigid_sphere_scattering(0.1, 1000.0, 343.0, 1.21);
529
530        assert!(problem.mesh.elements.len() > 0);
531        assert!(problem.mesh.nodes.nrows() > 0);
532        assert!(problem.ka() > 0.0);
533    }
534
535    #[test]
536    fn test_bem_solver_creation() {
537        let solver = BemSolver::new()
538            .with_solver_method(SolverMethod::Direct)
539            .with_assembly_method(AssemblyMethod::Tbem)
540            .with_verbose(false);
541
542        assert_eq!(solver.solver_method, SolverMethod::Direct);
543        assert_eq!(solver.assembly_method, AssemblyMethod::Tbem);
544    }
545
546    // Direct solver tests work in both native and WASM modes
547    #[test]
548    fn test_bem_solver_small_problem() {
549        // Very small problem for quick test
550        let problem = BemProblem::rigid_sphere_scattering_custom(
551            0.1,   // radius
552            100.0, // very low frequency for quick test
553            343.0, 1.21, 4, // coarse mesh
554            8,
555        );
556
557        let solver = BemSolver::new();
558        let result = solver.solve(&problem);
559
560        assert!(result.is_ok());
561        let solution = result.unwrap();
562        assert!(solution.num_dofs() > 0);
563        assert!(solution.max_surface_pressure() > 0.0);
564    }
565
566    #[test]
567    fn test_field_evaluation() {
568        // Very small problem
569        let problem = BemProblem::rigid_sphere_scattering_custom(0.1, 100.0, 343.0, 1.21, 4, 8);
570
571        let solver = BemSolver::new();
572        let solution = solver.solve(&problem).unwrap();
573
574        // Evaluate at a point outside the sphere
575        let p = solution.evaluate_pressure(&[0.0, 0.0, 0.2]);
576
577        // Should have some pressure (incident + scattered)
578        assert!(p.norm() > 0.0);
579    }
580}