1use 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#[derive(Debug, Clone, Copy, PartialEq, Eq, Default)]
48pub enum SolverMethod {
49 #[default]
51 Direct,
52 Cgs,
54 BiCgStab,
56}
57
58#[derive(Debug, Clone, Copy, PartialEq, Eq, Default)]
60pub enum AssemblyMethod {
61 #[default]
63 Tbem,
64 Slfmm,
66 Mlfmm,
68}
69
70#[derive(Debug, Clone, Copy, PartialEq, Eq)]
72pub enum BoundaryConditionType {
73 Rigid,
75 Soft,
77 Impedance,
79}
80
81#[derive(Debug, Clone)]
83pub struct BemProblem {
84 pub mesh: Mesh,
86 pub physics: PhysicsParams,
88 pub incident_field: IncidentField,
90 pub bc_type: BoundaryConditionType,
92 pub use_burton_miller: bool,
94}
95
96impl BemProblem {
97 pub fn rigid_sphere_scattering(
105 radius: f64,
106 frequency: f64,
107 speed_of_sound: f64,
108 density: f64,
109 ) -> Self {
110 let k = 2.0 * PI * frequency / speed_of_sound;
112 let ka = k * radius;
113
114 let subdivisions = if ka < 1.0 {
117 2 } else if ka < 5.0 {
119 3 } else {
121 4 };
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 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 pub fn with_incident_field(mut self, field: IncidentField) -> Self {
161 self.incident_field = field;
162 self
163 }
164
165 pub fn with_boundary_condition(mut self, bc_type: BoundaryConditionType) -> Self {
167 self.bc_type = bc_type;
168 self
169 }
170
171 pub fn with_burton_miller(mut self, use_bm: bool) -> Self {
173 self.use_burton_miller = use_bm;
174 self
175 }
176
177 pub fn ka(&self) -> f64 {
179 self.physics.wave_number * self.mesh_radius()
180 }
181
182 fn mesh_radius(&self) -> f64 {
184 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#[derive(Debug, Clone)]
199pub struct BemSolver {
200 pub solver_method: SolverMethod,
202 pub assembly_method: AssemblyMethod,
204 pub max_iterations: usize,
206 pub tolerance: f64,
208 pub verbose: bool,
210 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, }
224 }
225}
226
227impl BemSolver {
228 pub fn new() -> Self {
230 Self::default()
231 }
232
233 pub fn with_solver_method(mut self, method: SolverMethod) -> Self {
235 self.solver_method = method;
236 self
237 }
238
239 pub fn with_assembly_method(mut self, method: AssemblyMethod) -> Self {
241 self.assembly_method = method;
242 self
243 }
244
245 pub fn with_max_iterations(mut self, max_iter: usize) -> Self {
247 self.max_iterations = max_iter;
248 self
249 }
250
251 pub fn with_tolerance(mut self, tol: f64) -> Self {
253 self.tolerance = tol;
254 self
255 }
256
257 pub fn with_verbose(mut self, verbose: bool) -> Self {
259 self.verbose = verbose;
260 self
261 }
262
263 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 let elements = self.prepare_elements(problem);
281
282 let (matrix, rhs) =
284 self.assemble_system(&elements, &problem.mesh.nodes, &problem.physics)?;
285
286 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 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 fn prepare_elements(&self, problem: &BemProblem) -> Vec<Element> {
319 let mut elements = problem.mesh.elements.clone();
320
321 let bc = match problem.bc_type {
323 BoundaryConditionType::Rigid => {
324 BoundaryCondition::Velocity(vec![Complex64::new(0.0, 0.0)])
326 }
327 BoundaryConditionType::Soft => {
328 BoundaryCondition::Pressure(vec![Complex64::new(0.0, 0.0)])
330 }
331 BoundaryConditionType::Impedance => {
332 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 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 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 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 Err(BemError::NotImplemented(
367 "FMM assembly not yet available in high-level API".to_string(),
368 ))
369 }
370 }
371 }
372
373 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 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 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(¢ers, &normals, physics, beta)
398 } else {
399 incident_field.compute_rhs(¢ers, &normals, physics, false)
400 };
401
402 rhs = rhs + incident_rhs;
404
405 rhs
406 }
407
408 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 lu_solve(matrix, rhs).map_err(|e| BemError::SolverFailed(e.to_string()))
418 }
419 SolverMethod::Cgs | SolverMethod::BiCgStab => {
420 Err(BemError::NotImplemented(
422 "Iterative solvers not yet available in high-level API".to_string(),
423 ))
424 }
425 }
426 }
427}
428
429#[derive(Debug, Clone)]
431pub struct BemSolution {
432 pub surface_pressure: Array1<Complex64>,
434 pub elements: Vec<Element>,
436 pub nodes: Array2<f64>,
438 pub incident_field: IncidentField,
440 pub physics: PhysicsParams,
442}
443
444impl BemSolution {
445 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 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 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 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 pub fn num_dofs(&self) -> usize {
492 self.surface_pressure.len()
493 }
494}
495
496#[derive(Debug, Clone)]
498pub enum BemError {
499 NotImplemented(String),
501 SolverFailed(String),
503 InvalidMesh(String),
505 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 #[test]
548 fn test_bem_solver_small_problem() {
549 let problem = BemProblem::rigid_sphere_scattering_custom(
551 0.1, 100.0, 343.0, 1.21, 4, 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 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 let p = solution.evaluate_pressure(&[0.0, 0.0, 0.2]);
576
577 assert!(p.norm() > 0.0);
579 }
580}