1use super::*;
6use ndarray::{Array1, Array2};
7use num_complex::Complex64;
8use rayon::prelude::*;
9use std::f64::consts::PI;
10
11fn greens_function_3d(r: f64, k: f64) -> Complex64 {
14 if r < 1e-10 {
15 return Complex64::new(0.0, 0.0);
16 }
17 let ikr = Complex64::new(0.0, k * r);
18 ikr.exp() / (4.0 * PI * r)
19}
20
21fn greens_function_derivative(r: f64, k: f64, cos_angle: f64) -> Complex64 {
24 if r < 1e-10 {
25 return Complex64::new(0.0, 0.0);
26 }
27 let ikr = Complex64::new(0.0, k * r);
28 let factor = (ikr - 1.0) * ikr.exp() / (4.0 * PI * r * r);
29 factor * cos_angle
30}
31
32fn element_center_and_normal(nodes: &[Point3D]) -> (Point3D, Point3D) {
34 let center = Point3D::new(
36 nodes.iter().map(|n| n.x).sum::<f64>() / nodes.len() as f64,
37 nodes.iter().map(|n| n.y).sum::<f64>() / nodes.len() as f64,
38 nodes.iter().map(|n| n.z).sum::<f64>() / nodes.len() as f64,
39 );
40
41 let v1 = Point3D::new(
43 nodes[1].x - nodes[0].x,
44 nodes[1].y - nodes[0].y,
45 nodes[1].z - nodes[0].z,
46 );
47 let v2 = Point3D::new(
48 nodes[2].x - nodes[0].x,
49 nodes[2].y - nodes[0].y,
50 nodes[2].z - nodes[0].z,
51 );
52
53 let nx = v1.y * v2.z - v1.z * v2.y;
55 let ny = v1.z * v2.x - v1.x * v2.z;
56 let nz = v1.x * v2.y - v1.y * v2.x;
57
58 let norm = (nx * nx + ny * ny + nz * nz).sqrt();
59 let normal = Point3D::new(nx / norm, ny / norm, nz / norm);
60
61 (center, normal)
62}
63
64fn element_area(nodes: &[Point3D]) -> f64 {
66 if nodes.len() == 3 {
67 let v1 = Point3D::new(
69 nodes[1].x - nodes[0].x,
70 nodes[1].y - nodes[0].y,
71 nodes[1].z - nodes[0].z,
72 );
73 let v2 = Point3D::new(
74 nodes[2].x - nodes[0].x,
75 nodes[2].y - nodes[0].y,
76 nodes[2].z - nodes[0].z,
77 );
78
79 let cross_x = v1.y * v2.z - v1.z * v2.y;
80 let cross_y = v1.z * v2.x - v1.x * v2.z;
81 let cross_z = v1.x * v2.y - v1.y * v2.x;
82
83 0.5 * (cross_x * cross_x + cross_y * cross_y + cross_z * cross_z).sqrt()
84 } else if nodes.len() == 4 {
85 let v1 = Point3D::new(
87 nodes[1].x - nodes[0].x,
88 nodes[1].y - nodes[0].y,
89 nodes[1].z - nodes[0].z,
90 );
91 let v2 = Point3D::new(
92 nodes[2].x - nodes[0].x,
93 nodes[2].y - nodes[0].y,
94 nodes[2].z - nodes[0].z,
95 );
96
97 let cross1_x = v1.y * v2.z - v1.z * v2.y;
98 let cross1_y = v1.z * v2.x - v1.x * v2.z;
99 let cross1_z = v1.x * v2.y - v1.y * v2.x;
100 let area1 = 0.5 * (cross1_x * cross1_x + cross1_y * cross1_y + cross1_z * cross1_z).sqrt();
101
102 let v3 = Point3D::new(
103 nodes[3].x - nodes[0].x,
104 nodes[3].y - nodes[0].y,
105 nodes[3].z - nodes[0].z,
106 );
107
108 let cross2_x = v2.y * v3.z - v2.z * v3.y;
109 let cross2_y = v2.z * v3.x - v2.x * v3.z;
110 let cross2_z = v2.x * v3.y - v2.y * v3.x;
111 let area2 = 0.5 * (cross2_x * cross2_x + cross2_y * cross2_y + cross2_z * cross2_z).sqrt();
112
113 area1 + area2
114 } else {
115 0.0
116 }
117}
118
119pub fn build_bem_matrix(mesh: &RoomMesh, k: f64) -> Array2<Complex64> {
121 let n = mesh.elements.len();
122 let mut matrix = Array2::zeros((n, n));
123
124 let mut centers = Vec::new();
126 let mut normals = Vec::new();
127 let mut areas = Vec::new();
128
129 for element in &mesh.elements {
130 let nodes: Vec<Point3D> = element.nodes.iter().map(|&i| mesh.nodes[i]).collect();
131 let (center, normal) = element_center_and_normal(&nodes);
132 let area = element_area(&nodes);
133
134 centers.push(center);
135 normals.push(normal);
136 areas.push(area);
137 }
138
139 for i in 0..n {
142 for j in 0..n {
143 let r = centers[i].distance_to(¢ers[j]);
144
145 if i == j {
146 matrix[[i, j]] = Complex64::new(0.0, -k / (2.0 * PI)) * areas[j];
149 } else {
150 let dx = centers[i].x - centers[j].x;
152 let dy = centers[i].y - centers[j].y;
153 let dz = centers[i].z - centers[j].z;
154
155 let cos_angle = (dx * normals[i].x + dy * normals[i].y + dz * normals[i].z) / r;
157
158 matrix[[i, j]] = greens_function_derivative(r, k, cos_angle) * areas[j];
159 }
160 }
161 }
162
163 matrix
164}
165
166pub fn calculate_incident_field(
168 mesh: &RoomMesh,
169 sources: &[Source],
170 k: f64,
171 frequency: f64,
172) -> Array1<Complex64> {
173 let n = mesh.elements.len();
174 let mut incident = Array1::zeros(n);
175
176 for (i, element) in mesh.elements.iter().enumerate() {
177 let nodes: Vec<Point3D> = element.nodes.iter().map(|&idx| mesh.nodes[idx]).collect();
178 let (center, _normal) = element_center_and_normal(&nodes);
179
180 let mut total_pressure = Complex64::new(0.0, 0.0);
181
182 for source in sources {
183 let r = center.distance_to(&source.position);
184 let amplitude = source.amplitude_towards(¢er, frequency);
185
186 total_pressure += greens_function_3d(r, k) * amplitude;
188 }
189
190 incident[i] = total_pressure;
191 }
192
193 incident
194}
195
196pub fn calculate_field_pressure(
198 mesh: &RoomMesh,
199 surface_pressure: &Array1<Complex64>,
200 sources: &[Source],
201 field_points: &[Point3D],
202 k: f64,
203 frequency: f64,
204) -> Array1<Complex64> {
205 let n_points = field_points.len();
206 let mut pressures = Array1::zeros(n_points);
207
208 let mut centers = Vec::new();
210 let mut normals = Vec::new();
211 let mut areas = Vec::new();
212
213 for element in &mesh.elements {
214 let nodes: Vec<Point3D> = element.nodes.iter().map(|&i| mesh.nodes[i]).collect();
215 let (center, normal) = element_center_and_normal(&nodes);
216 let area = element_area(&nodes);
217
218 centers.push(center);
219 normals.push(normal);
220 areas.push(area);
221 }
222
223 for (ip, point) in field_points.iter().enumerate() {
224 let mut p_incident = Complex64::new(0.0, 0.0);
226 for source in sources {
227 let r = point.distance_to(&source.position);
228 let amplitude = source.amplitude_towards(point, frequency);
229 p_incident += greens_function_3d(r, k) * amplitude;
230 }
231
232 let mut p_scattered = Complex64::new(0.0, 0.0);
236
237 for j in 0..surface_pressure.len() {
238 let r = point.distance_to(¢ers[j]);
239 let g = greens_function_3d(r, k);
240
241 p_scattered += g * surface_pressure[j] * areas[j];
244 }
245
246 pressures[ip] = p_incident + p_scattered;
247 }
248
249 pressures
250}
251
252pub fn pressure_to_spl(pressure: Complex64) -> f64 {
254 let magnitude = pressure.norm();
255 let p_ref = 20e-6; 20.0 * (magnitude / p_ref).max(1e-10).log10()
257}
258
259pub fn gmres_solve(
262 a: &Array2<Complex64>,
263 b: &Array1<Complex64>,
264 max_iter: usize,
265 restart: usize,
266 tol: f64,
267) -> Result<Array1<Complex64>, String> {
268 let n = b.len();
269 if a.nrows() != n || a.ncols() != n {
270 return Err("Matrix dimensions mismatch".to_string());
271 }
272
273 let mut x = Array1::zeros(n);
274
275 for _cycle in 0..max_iter {
276 let ax = a.dot(&x);
278 let r = b - &ax;
279 let beta = r.iter().map(|ri| ri.norm_sqr()).sum::<f64>().sqrt();
280
281 if beta < tol {
282 return Ok(x);
283 }
284
285 let m = restart.min(n);
287 let mut v = vec![Array1::zeros(n); m + 1];
288 let mut h = Array2::<Complex64>::zeros((m + 1, m));
289
290 v[0] = r.mapv(|ri| ri / Complex64::new(beta, 0.0));
291
292 for j in 0..m {
293 let w = a.dot(&v[j]);
295 let mut w_orth = w.clone();
296
297 for i in 0..=j {
299 h[[i, j]] = v[i].iter().zip(w_orth.iter())
300 .map(|(vi, wi)| vi.conj() * wi)
301 .sum();
302
303 for k in 0..n {
304 w_orth[k] -= h[[i, j]] * v[i][k];
305 }
306 }
307
308 let h_norm = w_orth.iter().map(|wi| wi.norm_sqr()).sum::<f64>().sqrt();
309 h[[j + 1, j]] = Complex64::new(h_norm, 0.0);
310
311 if h_norm > 1e-12 {
312 v[j + 1] = w_orth.mapv(|wi| wi / Complex64::new(h_norm, 0.0));
313 } else {
314 break;
315 }
316 }
317
318 let mut e1 = Array1::<Complex64>::zeros(m + 1);
320 e1[0] = Complex64::new(beta, 0.0);
321
322 let y = solve_least_squares(&h, &e1, m)?;
324
325 for j in 0..m {
327 for k in 0..n {
328 x[k] += y[j] * v[j][k];
329 }
330 }
331
332 let ax = a.dot(&x);
334 let r_final = b - &ax;
335 let residual = r_final.iter().map(|ri| ri.norm_sqr()).sum::<f64>().sqrt();
336
337 if residual < tol {
338 return Ok(x);
339 }
340 }
341
342 Ok(x)
343}
344
345fn solve_least_squares(
347 h: &Array2<Complex64>,
348 e1: &Array1<Complex64>,
349 m: usize,
350) -> Result<Array1<Complex64>, String> {
351 let mut y = Array1::<Complex64>::zeros(m);
352 let mut rhs = e1.slice(ndarray::s![0..m]).to_owned();
353
354 let mut h_tri = h.slice(ndarray::s![0..m, 0..m]).to_owned();
356
357 for i in 0..m {
358 for j in (i + 1)..m {
359 if h_tri[[j, i]].norm() > 1e-12 {
360 let a = h_tri[[i, i]];
361 let b = h_tri[[j, i]];
362 let r = (a.norm_sqr() + b.norm_sqr()).sqrt();
363 let c = a.norm() / r;
364 let s = b / Complex64::new(r, 0.0);
365
366 for k in i..m {
368 let temp = c * h_tri[[i, k]] + s * h_tri[[j, k]];
369 h_tri[[j, k]] = -s.conj() * h_tri[[i, k]] + c * h_tri[[j, k]];
370 h_tri[[i, k]] = temp;
371 }
372
373 let temp = c * rhs[i] + s * rhs[j];
374 rhs[j] = -s.conj() * rhs[i] + c * rhs[j];
375 rhs[i] = temp;
376 }
377 }
378 }
379
380 for i in (0..m).rev() {
382 let mut sum = rhs[i];
383 for j in (i + 1)..m {
384 sum -= h_tri[[i, j]] * y[j];
385 }
386 if h_tri[[i, i]].norm() > 1e-12 {
387 y[i] = sum / h_tri[[i, i]];
388 }
389 }
390
391 Ok(y)
392}
393
394pub fn solve_bem_system(
396 mesh: &RoomMesh,
397 sources: &[Source],
398 k: f64,
399 frequency: f64,
400) -> Result<Array1<Complex64>, String> {
401 println!(" Building BEM matrix (parallel)...");
402 let matrix = build_bem_matrix_parallel(mesh, k);
403
404 println!(" Computing RHS (parallel)...");
405 let rhs = calculate_incident_field_derivative_parallel(mesh, sources, k, frequency);
406
407 println!(" Solving with GMRES...");
408 let max_iter = 100;
409 let restart = 50;
410 let tol = 1e-6;
411
412 gmres_solve(&matrix, &rhs, max_iter, restart, tol)
413}
414
415pub fn build_bem_matrix_parallel(mesh: &RoomMesh, k: f64) -> Array2<Complex64> {
417 let n = mesh.elements.len();
418
419 let element_data: Vec<_> = mesh.elements.par_iter()
421 .map(|element| {
422 let nodes: Vec<Point3D> = element.nodes.iter().map(|&i| mesh.nodes[i]).collect();
423 let (center, normal) = element_center_and_normal(&nodes);
424 let area = element_area(&nodes);
425 (center, normal, area)
426 })
427 .collect();
428
429 let rows: Vec<_> = (0..n).into_par_iter()
431 .map(|i| {
432 let mut row = vec![Complex64::new(0.0, 0.0); n];
433 let (center_i, normal_i, _area_i) = &element_data[i];
434
435 for j in 0..n {
436 let (center_j, _normal_j, area_j) = &element_data[j];
437 let r = center_i.distance_to(center_j);
438
439 if i == j {
440 row[j] = Complex64::new(0.0, -k / (2.0 * PI)) * area_j;
442 } else {
443 let dx = center_i.x - center_j.x;
445 let dy = center_i.y - center_j.y;
446 let dz = center_i.z - center_j.z;
447
448 let cos_angle = (dx * normal_i.x + dy * normal_i.y + dz * normal_i.z) / r;
449 row[j] = greens_function_derivative(r, k, cos_angle) * area_j;
450 }
451 }
452 row
453 })
454 .collect();
455
456 let mut matrix = Array2::zeros((n, n));
458 for (i, row) in rows.iter().enumerate() {
459 for (j, &val) in row.iter().enumerate() {
460 matrix[[i, j]] = val;
461 }
462 }
463
464 matrix
465}
466
467pub fn build_bem_matrix_adaptive(mesh: &RoomMesh, k: f64, use_adaptive: bool) -> Array2<Complex64> {
473 use crate::core::integration::singular::QuadratureParams;
474 use crate::core::types::{ElementType, PhysicsParams};
475
476 let n = mesh.elements.len();
477
478 let speed_of_sound = 343.0;
480 let frequency = k * speed_of_sound / (2.0 * PI);
481 let omega = 2.0 * PI * frequency;
482
483 let physics = PhysicsParams {
484 wave_number: k,
485 density: 1.0,
486 speed_of_sound,
487 frequency,
488 omega,
489 wave_length: speed_of_sound / frequency,
490 harmonic_factor: 1.0, pressure_factor: 1.0 * omega * 1.0,
492 tau: -1.0, };
494
495 let element_data: Vec<_> = mesh.elements.par_iter()
497 .map(|element| {
498 let nodes: Vec<Point3D> = element.nodes.iter().map(|&i| mesh.nodes[i]).collect();
499 let (center, normal) = element_center_and_normal(&nodes);
500 let area = element_area(&nodes);
501 let char_length = element_characteristic_length(&nodes);
502 (center, normal, area, char_length, nodes)
503 })
504 .collect();
505
506 let rows: Vec<_> = (0..n).into_par_iter()
508 .map(|i| {
509 let mut row = vec![Complex64::new(0.0, 0.0); n];
510 let (center_i, normal_i, _area_i, char_length_i, _nodes_i) = &element_data[i];
511
512 for j in 0..n {
513 let (center_j, _normal_j, area_j, char_length_j, nodes_j) = &element_data[j];
514 let r = center_i.distance_to(center_j);
515
516 let near_threshold = 2.0 * (char_length_i + char_length_j);
518 let is_near = r < near_threshold || i == j;
519
520 if use_adaptive && is_near {
521 let element_coords = nodes_to_array2(nodes_j);
523 let source_point = point_to_array1(center_i);
524 let source_normal = normal_to_array1(normal_i);
525
526 let ka = k * char_length_j;
528 let quad_params = QuadratureParams::for_ka(ka);
529
530 let result = crate::core::integration::singular::singular_integration_with_params(
531 &source_point,
532 &source_normal,
533 &element_coords,
534 ElementType::Tri3,
535 &physics,
536 None,
537 0, false, &quad_params,
540 );
541
542 row[j] = result.dg_dn_integral;
544 } else {
545 if i == j {
547 row[j] = Complex64::new(0.0, -k / (2.0 * PI)) * area_j;
549 } else {
550 let dx = center_i.x - center_j.x;
552 let dy = center_i.y - center_j.y;
553 let dz = center_i.z - center_j.z;
554
555 let cos_angle = (dx * normal_i.x + dy * normal_i.y + dz * normal_i.z) / r;
556 row[j] = greens_function_derivative(r, k, cos_angle) * area_j;
557 }
558 }
559 }
560 row
561 })
562 .collect();
563
564 let mut matrix = Array2::zeros((n, n));
566 for (i, row) in rows.iter().enumerate() {
567 for (j, &val) in row.iter().enumerate() {
568 matrix[[i, j]] = val;
569 }
570 }
571
572 matrix
573}
574
575fn element_characteristic_length(nodes: &[Point3D]) -> f64 {
577 if nodes.len() < 3 {
578 return 0.0;
579 }
580
581 let d01 = nodes[0].distance_to(&nodes[1]);
583 let d12 = nodes[1].distance_to(&nodes[2]);
584 let d20 = nodes[2].distance_to(&nodes[0]);
585
586 (d01 + d12 + d20) / 3.0
587}
588
589fn point_to_array1(p: &Point3D) -> Array1<f64> {
591 use ndarray::array;
592 array![p.x, p.y, p.z]
593}
594
595fn normal_to_array1(n: &Point3D) -> Array1<f64> {
597 use ndarray::array;
598 array![n.x, n.y, n.z]
599}
600
601fn nodes_to_array2(nodes: &[Point3D]) -> Array2<f64> {
603 let n = nodes.len();
604 let mut coords = Array2::zeros((n, 3));
605 for (i, node) in nodes.iter().enumerate() {
606 coords[[i, 0]] = node.x;
607 coords[[i, 1]] = node.y;
608 coords[[i, 2]] = node.z;
609 }
610 coords
611}
612
613pub fn calculate_incident_field_derivative_parallel(
615 mesh: &RoomMesh,
616 sources: &[Source],
617 k: f64,
618 frequency: f64,
619) -> Array1<Complex64> {
620 let element_data: Vec<_> = mesh.elements.par_iter()
621 .map(|element| {
622 let nodes: Vec<Point3D> = element.nodes.iter().map(|&idx| mesh.nodes[idx]).collect();
623 let (center, normal) = element_center_and_normal(&nodes);
624
625 let mut dpdn_inc = Complex64::new(0.0, 0.0);
627
628 for source in sources {
629 let r = center.distance_to(&source.position);
630 if r < 1e-10 {
631 continue;
632 }
633
634 let amplitude = source.amplitude_towards(¢er, frequency);
635
636 let dx = center.x - source.position.x;
638 let dy = center.y - source.position.y;
639 let dz = center.z - source.position.z;
640
641 let cos_angle = (dx * normal.x + dy * normal.y + dz * normal.z) / r;
643
644 dpdn_inc += greens_function_derivative(r, k, cos_angle) * amplitude;
645 }
646
647 -dpdn_inc
652 })
653 .collect();
654
655 Array1::from_vec(element_data)
656}
657
658pub fn calculate_field_pressure_bem_parallel(
660 mesh: &RoomMesh,
661 surface_normal_derivative: &Array1<Complex64>,
662 sources: &[Source],
663 field_points: &[Point3D],
664 k: f64,
665 frequency: f64,
666) -> Array1<Complex64> {
667 let element_data: Vec<_> = mesh.elements.iter()
669 .map(|element| {
670 let nodes: Vec<Point3D> = element.nodes.iter().map(|&i| mesh.nodes[i]).collect();
671 let (center, _normal) = element_center_and_normal(&nodes);
672 let area = element_area(&nodes);
673 (center, area)
674 })
675 .collect();
676
677 let pressures: Vec<_> = field_points.par_iter()
679 .map(|point| {
680 let mut p_incident = Complex64::new(0.0, 0.0);
682 for source in sources {
683 let r = point.distance_to(&source.position);
684 if r < 1e-10 {
685 continue;
686 }
687 let amplitude = source.amplitude_towards(point, frequency);
688 p_incident += greens_function_3d(r, k) * amplitude;
689 }
690
691 let mut p_scattered = Complex64::new(0.0, 0.0);
694 for (j, (center_j, area_j)) in element_data.iter().enumerate() {
695 let r = point.distance_to(center_j);
696 if r < 1e-10 {
697 continue;
698 }
699 let g = greens_function_3d(r, k);
700 p_scattered += g * surface_normal_derivative[j] * area_j;
701 }
702
703 p_incident + p_scattered
704 })
705 .collect();
706
707 Array1::from_vec(pressures)
708}
709
710use crate::core::mesh::octree::Octree;
715use crate::core::types::{
716 BoundaryCondition, Cluster, Element, ElementProperty, ElementType, PhysicsParams,
717};
718use crate::core::assembly::slfmm::{build_slfmm_system, SlfmmSystem};
719use crate::core::solver::{
720 gmres_solve_with_ilu_operator, GmresConfig, GmresSolution,
721 IluMethod, IluScanningDegree, SlfmmOperator,
722};
723
724pub struct FmmSolverConfig {
726 pub max_elements_per_leaf: usize,
728 pub max_tree_depth: usize,
730 pub n_theta: usize,
732 pub n_phi: usize,
734 pub n_terms: usize,
736 pub separation_ratio: f64,
738}
739
740impl Default for FmmSolverConfig {
741 fn default() -> Self {
742 Self {
743 max_elements_per_leaf: 50,
744 max_tree_depth: 8,
745 n_theta: 6,
746 n_phi: 12,
747 n_terms: 6,
748 separation_ratio: 1.5, }
750 }
751}
752
753pub fn room_mesh_to_core_elements(
755 mesh: &RoomMesh,
756 _k: f64,
757) -> (Vec<Element>, Array2<f64>) {
758 let n_nodes = mesh.nodes.len();
759 let n_elements = mesh.elements.len();
760
761 let mut nodes = Array2::zeros((n_nodes, 3));
763 for (i, node) in mesh.nodes.iter().enumerate() {
764 nodes[[i, 0]] = node.x;
765 nodes[[i, 1]] = node.y;
766 nodes[[i, 2]] = node.z;
767 }
768
769 let mut elements = Vec::with_capacity(n_elements);
771 for (elem_idx, surface_elem) in mesh.elements.iter().enumerate() {
772 let elem_nodes: Vec<Point3D> = surface_elem.nodes.iter()
773 .map(|&i| mesh.nodes[i])
774 .collect();
775
776 let (center, normal) = element_center_and_normal(&elem_nodes);
777 let area = element_area(&elem_nodes);
778
779 let elem_type = if surface_elem.nodes.len() == 3 {
781 ElementType::Tri3
782 } else {
783 ElementType::Quad4
784 };
785
786 let element = Element {
787 connectivity: surface_elem.nodes.clone(),
788 element_type: elem_type,
789 property: ElementProperty::Surface,
790 normal: Array1::from_vec(vec![normal.x, normal.y, normal.z]),
791 node_normals: Array2::zeros((surface_elem.nodes.len(), 3)),
792 center: Array1::from_vec(vec![center.x, center.y, center.z]),
793 area,
794 boundary_condition: BoundaryCondition::Velocity(vec![Complex64::new(0.0, 0.0)]),
795 group: 0,
796 dof_addresses: vec![elem_idx],
797 };
798
799 elements.push(element);
800 }
801
802 (elements, nodes)
803}
804
805pub fn build_clusters_from_octree(
807 octree: &Octree,
808 elements: &[Element],
809) -> Vec<Cluster> {
810 let leaves = octree.leaves();
811 let mut clusters = Vec::with_capacity(leaves.len());
812
813 for &leaf_idx in &leaves {
814 let node = &octree.nodes[leaf_idx];
815 if node.element_indices.is_empty() {
816 continue;
817 }
818
819 let mut cluster = Cluster::new(node.center.clone());
820 cluster.element_indices = node.element_indices.clone();
821 cluster.num_elements = node.element_indices.len();
822 cluster.element_property = ElementProperty::Surface;
823 cluster.radius = node.radius();
824 cluster.level = node.level;
825
826 cluster.num_dofs = node.element_indices.iter()
828 .filter(|&&i| !elements[i].property.is_evaluation())
829 .count();
830 cluster.dofs_per_element = 1;
831
832 clusters.push(cluster);
833 }
834
835 let leaf_to_cluster: std::collections::HashMap<usize, usize> = leaves.iter()
838 .filter(|&&i| !octree.nodes[i].element_indices.is_empty())
839 .enumerate()
840 .map(|(cluster_idx, &leaf_idx)| (leaf_idx, cluster_idx))
841 .collect();
842
843 for (cluster_idx, &leaf_idx) in leaves.iter().enumerate() {
845 if octree.nodes[leaf_idx].element_indices.is_empty() {
846 continue;
847 }
848
849 let octree_node = &octree.nodes[leaf_idx];
850
851 let near: Vec<usize> = octree_node.near_clusters.iter()
853 .filter_map(|&near_leaf| leaf_to_cluster.get(&near_leaf).copied())
854 .collect();
855
856 let far: Vec<usize> = octree_node.far_clusters.iter()
858 .filter_map(|&far_leaf| leaf_to_cluster.get(&far_leaf).copied())
859 .collect();
860
861 if let Some(cluster) = clusters.get_mut(cluster_idx) {
862 cluster.near_clusters = near;
863 cluster.far_clusters = far;
864 }
865 }
866
867 clusters
868}
869
870pub fn build_fmm_system(
872 mesh: &RoomMesh,
873 sources: &[Source],
874 k: f64,
875 frequency: f64,
876 fmm_config: &FmmSolverConfig,
877) -> Result<(SlfmmSystem, Vec<Element>, Array2<f64>), String> {
878 println!(" Converting room mesh to core elements...");
879 let (elements, nodes) = room_mesh_to_core_elements(mesh, k);
880 println!(" {} elements, {} nodes", elements.len(), nodes.nrows());
881
882 println!(" Building octree...");
884 let centers: Vec<Array1<f64>> = elements.iter()
885 .map(|e| e.center.clone())
886 .collect();
887
888 let mut octree = Octree::build(
889 ¢ers,
890 fmm_config.max_elements_per_leaf,
891 fmm_config.max_tree_depth,
892 );
893
894 octree.compute_interaction_lists(fmm_config.separation_ratio);
896
897 let stats = octree.stats();
898 println!(" {} leaves, {} levels, avg {:.1} elements/leaf",
899 stats.num_leaves, stats.num_levels, stats.avg_elements_per_leaf);
900
901 println!(" Building clusters...");
903 let clusters = build_clusters_from_octree(&octree, &elements);
904 println!(" {} clusters", clusters.len());
905
906 let speed_of_sound = 343.0;
908 let physics = PhysicsParams::new(frequency, speed_of_sound, 1.21, true);
909
910 println!(" Assembling SLFMM system...");
912 let mut system = build_slfmm_system(
913 &elements,
914 &nodes,
915 &clusters,
916 &physics,
917 fmm_config.n_theta,
918 fmm_config.n_phi,
919 fmm_config.n_terms,
920 );
921
922 println!(" Computing RHS...");
924 let rhs = calculate_incident_field_derivative_parallel(mesh, sources, k, frequency);
925 system.rhs = rhs;
926
927 Ok((system, elements, nodes))
928}
929
930pub fn solve_bem_fmm_gmres_ilu(
938 mesh: &RoomMesh,
939 sources: &[Source],
940 k: f64,
941 frequency: f64,
942 fmm_config: &FmmSolverConfig,
943 gmres_max_iter: usize,
944 gmres_restart: usize,
945 gmres_tolerance: f64,
946 ilu_method: IluMethod,
947 ilu_degree: IluScanningDegree,
948) -> Result<Array1<Complex64>, String> {
949 let (system, _elements, _nodes) = build_fmm_system(mesh, sources, k, frequency, fmm_config)?;
951
952 println!(" Extracting near-field matrix for ILU...");
955 let nearfield_matrix = system.extract_near_field_matrix();
956 println!(" Near-field matrix: {}x{}", nearfield_matrix.nrows(), nearfield_matrix.ncols());
957
958 let fmm_operator = SlfmmOperator::new(system);
960
961 let gmres_config = GmresConfig {
963 max_iterations: gmres_max_iter,
964 restart: gmres_restart,
965 tolerance: gmres_tolerance,
966 print_interval: 0,
967 };
968
969 let rhs = fmm_operator.rhs().clone();
971
972 println!(" Solving with FMM + GMRES + ILU...");
974 let result = gmres_solve_with_ilu_operator(
975 &fmm_operator,
976 &nearfield_matrix,
977 &rhs,
978 ilu_method,
979 ilu_degree,
980 &gmres_config,
981 );
982
983 if result.converged {
984 println!(" Converged in {} iterations, residual: {:.2e}",
985 result.iterations, result.residual);
986 } else {
987 println!(" Warning: Did not converge after {} iterations, residual: {:.2e}",
988 result.iterations, result.residual);
989 }
990
991 Ok(result.x)
992}
993
994pub fn solve_bem_fmm_gmres_hierarchical(
1008 mesh: &RoomMesh,
1009 sources: &[Source],
1010 k: f64,
1011 frequency: f64,
1012 fmm_config: &FmmSolverConfig,
1013 gmres_max_iter: usize,
1014 gmres_restart: usize,
1015 gmres_tolerance: f64,
1016) -> Result<Array1<Complex64>, String> {
1017 use crate::core::solver::gmres_solve_with_hierarchical_precond;
1018
1019 let (system, _elements, _nodes) = build_fmm_system(mesh, sources, k, frequency, fmm_config)?;
1021
1022 let gmres_config = GmresConfig {
1024 max_iterations: gmres_max_iter,
1025 restart: gmres_restart,
1026 tolerance: gmres_tolerance,
1027 print_interval: 0,
1028 };
1029
1030 let rhs = system.rhs.clone();
1032
1033 println!(" Solving with FMM + GMRES + Hierarchical Preconditioner...");
1035 let result = gmres_solve_with_hierarchical_precond(&system, &rhs, &gmres_config);
1036
1037 if result.converged {
1038 println!(" Converged in {} iterations, residual: {:.2e}",
1039 result.iterations, result.residual);
1040 } else {
1041 println!(" Warning: Did not converge after {} iterations, residual: {:.2e}",
1042 result.iterations, result.residual);
1043 }
1044
1045 Ok(result.x)
1046}
1047
1048pub fn solve_bem_fmm_gmres_ilu_with_result(
1053 mesh: &RoomMesh,
1054 sources: &[Source],
1055 k: f64,
1056 frequency: f64,
1057 fmm_config: &FmmSolverConfig,
1058 gmres_max_iter: usize,
1059 gmres_restart: usize,
1060 gmres_tolerance: f64,
1061 ilu_method: IluMethod,
1062 ilu_degree: IluScanningDegree,
1063) -> Result<GmresSolution, String> {
1064 let (system, _elements, _nodes) = build_fmm_system(mesh, sources, k, frequency, fmm_config)?;
1066
1067 println!(" Extracting near-field matrix for ILU...");
1069 let nearfield_matrix = system.extract_near_field_matrix();
1070 println!(" Near-field matrix: {}x{}", nearfield_matrix.nrows(), nearfield_matrix.ncols());
1071
1072 let fmm_operator = SlfmmOperator::new(system);
1074
1075 let gmres_config = GmresConfig {
1077 max_iterations: gmres_max_iter,
1078 restart: gmres_restart,
1079 tolerance: gmres_tolerance,
1080 print_interval: 0,
1081 };
1082
1083 let rhs = fmm_operator.rhs().clone();
1085
1086 println!(" Solving with FMM + GMRES + ILU...");
1088 let result = gmres_solve_with_ilu_operator(
1089 &fmm_operator,
1090 &nearfield_matrix,
1091 &rhs,
1092 ilu_method,
1093 ilu_degree,
1094 &gmres_config,
1095 );
1096
1097 if result.converged {
1098 println!(" Converged in {} iterations, residual: {:.2e}",
1099 result.iterations, result.residual);
1100 } else {
1101 println!(" Warning: Did not converge after {} iterations, residual: {:.2e}",
1102 result.iterations, result.residual);
1103 }
1104
1105 Ok(result)
1106}
1107
1108#[cfg(test)]
1109mod tests {
1110 use super::*;
1111
1112 #[test]
1113 fn test_greens_function() {
1114 let k = 2.0 * PI * 1000.0 / 343.0;
1115 let r = 1.0;
1116 let g = greens_function_3d(r, k);
1117 assert!((g.norm() - 1.0 / (4.0 * PI)).abs() < 0.1);
1119 }
1120
1121 #[test]
1122 fn test_pressure_to_spl() {
1123 let p = Complex64::new(1.0, 0.0); let spl = pressure_to_spl(p);
1125 assert!((spl - 94.0).abs() < 1.0);
1127 }
1128
1129 #[test]
1130 fn test_room_mesh_to_core_elements() {
1131 let nodes = vec![
1133 Point3D::new(0.0, 0.0, 0.0),
1134 Point3D::new(1.0, 0.0, 0.0),
1135 Point3D::new(0.5, 1.0, 0.0),
1136 ];
1137 let elements = vec![
1138 SurfaceElement { nodes: vec![0, 1, 2] },
1139 ];
1140 let mesh = RoomMesh { nodes, elements };
1141
1142 let k = 2.0 * PI * 100.0 / 343.0;
1143 let (core_elements, core_nodes) = room_mesh_to_core_elements(&mesh, k);
1144
1145 assert_eq!(core_elements.len(), 1);
1146 assert_eq!(core_nodes.nrows(), 3);
1147 assert_eq!(core_elements[0].element_type, ElementType::Tri3);
1148 }
1149}