1use ndarray::{Array1, Array2};
21use num_complex::Complex64;
22
23use crate::core::greens::legendre::legendre_polynomials;
24use crate::core::greens::spherical::spherical_hankel_first_kind;
25use crate::core::integration::{regular_integration, singular_integration, unit_sphere_quadrature};
26use crate::core::types::{Cluster, ClusterLevel, Element, PhysicsParams};
27
28pub struct MlfmmSystem {
30 pub num_levels: usize,
32 pub near_matrix: Vec<NearFieldBlock>,
34 pub t_matrices: Vec<LevelMatrices>,
37 pub d_matrices: Vec<Vec<DMatrixEntry>>,
39 pub s_matrices: Vec<LevelMatrices>,
41 pub rhs: Array1<Complex64>,
43 pub num_dofs: usize,
45 cluster_levels: Vec<ClusterLevel>,
47 pub leaf_dof_indices: Vec<Vec<usize>>,
49 pub sphere_points_per_level: Vec<usize>,
51}
52
53#[derive(Debug, Clone)]
55pub struct NearFieldBlock {
56 pub source_cluster: usize,
58 pub field_cluster: usize,
60 pub coefficients: Array2<Complex64>,
62}
63
64#[derive(Debug, Clone)]
66pub struct DMatrixEntry {
67 pub source_cluster: usize,
69 pub field_cluster: usize,
71 pub level: usize,
73 pub coefficients: Array2<Complex64>,
75}
76
77#[derive(Debug, Clone)]
79pub struct LevelMatrices {
80 pub level: usize,
82 pub num_clusters: usize,
84 pub matrices: Vec<Array2<Complex64>>,
86}
87
88impl MlfmmSystem {
89 pub fn new(num_dofs: usize, num_levels: usize, cluster_levels: Vec<ClusterLevel>) -> Self {
91 let sphere_points_per_level: Vec<usize> = cluster_levels
92 .iter()
93 .map(|level| level.theta_points * level.phi_points)
94 .collect();
95
96 Self {
97 num_levels,
98 near_matrix: Vec::new(),
99 t_matrices: Vec::with_capacity(num_levels),
100 d_matrices: Vec::with_capacity(num_levels),
101 s_matrices: Vec::with_capacity(num_levels),
102 rhs: Array1::zeros(num_dofs),
103 num_dofs,
104 cluster_levels,
105 leaf_dof_indices: Vec::new(),
106 sphere_points_per_level,
107 }
108 }
109
110 pub fn matvec(&self, x: &Array1<Complex64>) -> Array1<Complex64> {
127 let mut y = Array1::zeros(self.num_dofs);
128
129 if self.num_levels == 0 {
130 return y;
131 }
132
133 let leaf_level = self.num_levels - 1;
134
135 for block in &self.near_matrix {
141 if block.source_cluster >= self.leaf_dof_indices.len()
142 || block.field_cluster >= self.leaf_dof_indices.len()
143 {
144 continue;
145 }
146
147 let src_dofs = &self.leaf_dof_indices[block.source_cluster];
148 let fld_dofs = &self.leaf_dof_indices[block.field_cluster];
149
150 let x_local: Array1<Complex64> =
152 Array1::from_iter(fld_dofs.iter().map(|&i| x[i]));
153
154 let y_local = block.coefficients.dot(&x_local);
156
157 for (local_i, &global_i) in src_dofs.iter().enumerate() {
159 if local_i < y_local.len() {
160 y[global_i] += y_local[local_i];
161 }
162 }
163
164 if block.source_cluster != block.field_cluster {
167 let x_src: Array1<Complex64> =
169 Array1::from_iter(src_dofs.iter().map(|&i| x[i]));
170 let y_fld = block.coefficients.t().dot(&x_src);
172 for (local_j, &global_j) in fld_dofs.iter().enumerate() {
174 if local_j < y_fld.len() {
175 y[global_j] += y_fld[local_j];
176 }
177 }
178 }
179 }
180
181 if self.num_levels > 1 && !self.t_matrices.is_empty() {
183 let multipoles = self.upward_pass(x);
185
186 let mut locals = self.translate_all_levels(&multipoles);
188
189 self.downward_pass(&mut locals);
191
192 self.evaluate_locals(&locals, leaf_level, &mut y);
194 }
195
196 y
197 }
198
199 fn upward_pass(&self, x: &Array1<Complex64>) -> Vec<Vec<Array1<Complex64>>> {
204 let mut multipoles: Vec<Vec<Array1<Complex64>>> = Vec::with_capacity(self.num_levels);
205
206 for level in 0..self.num_levels {
208 let num_clusters = if level < self.cluster_levels.len() {
209 self.cluster_levels[level].clusters.len()
210 } else {
211 0
212 };
213 let num_points = if level < self.sphere_points_per_level.len() {
214 self.sphere_points_per_level[level]
215 } else {
216 0
217 };
218 multipoles.push(
219 (0..num_clusters)
220 .map(|_| Array1::zeros(num_points))
221 .collect(),
222 );
223 }
224
225 let leaf_level = self.num_levels - 1;
226
227 if leaf_level < self.t_matrices.len() {
229 let t_level = &self.t_matrices[leaf_level];
230 for cluster_idx in 0..t_level.num_clusters {
231 if cluster_idx >= self.leaf_dof_indices.len() {
232 continue;
233 }
234 let cluster_dofs = &self.leaf_dof_indices[cluster_idx];
235 if cluster_dofs.is_empty() {
236 continue;
237 }
238 let t_mat = &t_level.matrices[cluster_idx];
239 if t_mat.is_empty() {
240 continue;
241 }
242 let x_local: Array1<Complex64> =
243 Array1::from_iter(cluster_dofs.iter().map(|&i| x[i]));
244 if x_local.len() == t_mat.ncols() {
245 multipoles[leaf_level][cluster_idx] = t_mat.dot(&x_local);
246 }
247 }
248 }
249
250 for level in (0..leaf_level).rev() {
253 if level >= self.t_matrices.len() || level >= self.cluster_levels.len() {
254 continue;
255 }
256
257 let t_level = &self.t_matrices[level];
258 let clusters = &self.cluster_levels[level].clusters;
259
260 for (cluster_idx, cluster) in clusters.iter().enumerate() {
261 if cluster_idx >= t_level.matrices.len() {
262 continue;
263 }
264
265 let t_mat = &t_level.matrices[cluster_idx];
267 if t_mat.is_empty() || cluster.sons.is_empty() {
268 continue;
269 }
270
271 let child_level = level + 1;
273 if child_level >= multipoles.len() {
274 continue;
275 }
276
277 let child_num_points = self.sphere_points_per_level.get(child_level).copied().unwrap_or(0);
278 let mut child_multipoles = Array1::zeros(cluster.sons.len() * child_num_points);
279
280 for (son_idx, &child_cluster_idx) in cluster.sons.iter().enumerate() {
281 if child_cluster_idx < multipoles[child_level].len() {
282 let child_mult = &multipoles[child_level][child_cluster_idx];
283 let offset = son_idx * child_num_points;
284 for (i, &val) in child_mult.iter().enumerate() {
285 if offset + i < child_multipoles.len() {
286 child_multipoles[offset + i] = val;
287 }
288 }
289 }
290 }
291
292 if child_multipoles.len() == t_mat.ncols() {
293 multipoles[level][cluster_idx] = t_mat.dot(&child_multipoles);
294 }
295 }
296 }
297
298 multipoles
299 }
300
301 fn translate_all_levels(
303 &self,
304 multipoles: &[Vec<Array1<Complex64>>],
305 ) -> Vec<Vec<Array1<Complex64>>> {
306 let mut locals: Vec<Vec<Array1<Complex64>>> = Vec::with_capacity(self.num_levels);
307
308 for level in 0..self.num_levels {
310 let num_clusters = if level < self.cluster_levels.len() {
311 self.cluster_levels[level].clusters.len()
312 } else {
313 0
314 };
315 let num_points = self.sphere_points_per_level.get(level).copied().unwrap_or(0);
316 locals.push(
317 (0..num_clusters)
318 .map(|_| Array1::zeros(num_points))
319 .collect(),
320 );
321 }
322
323 for (level, d_entries) in self.d_matrices.iter().enumerate() {
325 if level >= multipoles.len() || level >= locals.len() {
326 continue;
327 }
328
329 for d_entry in d_entries {
330 if d_entry.source_cluster >= multipoles[level].len()
331 || d_entry.field_cluster >= locals[level].len()
332 {
333 continue;
334 }
335
336 let src_mult = &multipoles[level][d_entry.source_cluster];
337 if src_mult.len() != d_entry.coefficients.ncols() {
338 continue;
339 }
340
341 let translated = d_entry.coefficients.dot(src_mult);
342 for (i, &val) in translated.iter().enumerate() {
343 if i < locals[level][d_entry.field_cluster].len() {
344 locals[level][d_entry.field_cluster][i] += val;
345 }
346 }
347 }
348 }
349
350 locals
351 }
352
353 fn downward_pass(&self, locals: &mut [Vec<Array1<Complex64>>]) {
355 for level in 0..self.num_levels.saturating_sub(1) {
357 if level >= self.s_matrices.len() || level >= self.cluster_levels.len() {
358 continue;
359 }
360
361 let s_level = &self.s_matrices[level];
362 let clusters = &self.cluster_levels[level].clusters;
363
364 for (cluster_idx, cluster) in clusters.iter().enumerate() {
365 if cluster_idx >= s_level.matrices.len() || cluster.sons.is_empty() {
366 continue;
367 }
368
369 let s_mat = &s_level.matrices[cluster_idx];
370 if s_mat.is_empty() {
371 continue;
372 }
373
374 let parent_local = &locals[level][cluster_idx];
376 if parent_local.len() != s_mat.ncols() {
377 continue;
378 }
379
380 let child_locals = s_mat.dot(parent_local);
381 let child_level = level + 1;
382 if child_level >= locals.len() {
383 continue;
384 }
385
386 let child_num_points = self.sphere_points_per_level.get(child_level).copied().unwrap_or(0);
387
388 for (son_idx, &child_cluster_idx) in cluster.sons.iter().enumerate() {
390 if child_cluster_idx >= locals[child_level].len() {
391 continue;
392 }
393
394 let offset = son_idx * child_num_points;
395 for i in 0..child_num_points {
396 if offset + i < child_locals.len() {
397 if i < locals[child_level][child_cluster_idx].len() {
398 locals[child_level][child_cluster_idx][i] += child_locals[offset + i];
399 }
400 }
401 }
402 }
403 }
404 }
405 }
406
407 fn evaluate_locals(
409 &self,
410 locals: &[Vec<Array1<Complex64>>],
411 leaf_level: usize,
412 y: &mut Array1<Complex64>,
413 ) {
414 if leaf_level >= self.s_matrices.len() || leaf_level >= locals.len() {
415 return;
416 }
417
418 let s_level = &self.s_matrices[leaf_level];
419
420 for cluster_idx in 0..s_level.num_clusters {
421 if cluster_idx >= self.leaf_dof_indices.len() || cluster_idx >= locals[leaf_level].len() {
422 continue;
423 }
424
425 let cluster_dofs = &self.leaf_dof_indices[cluster_idx];
426 if cluster_dofs.is_empty() {
427 continue;
428 }
429
430 let s_mat = &s_level.matrices[cluster_idx];
431 if s_mat.is_empty() {
432 continue;
433 }
434
435 let local_exp = &locals[leaf_level][cluster_idx];
436 if local_exp.len() != s_mat.ncols() {
437 continue;
438 }
439
440 let y_local = s_mat.dot(local_exp);
441 for (local_j, &global_j) in cluster_dofs.iter().enumerate() {
442 if local_j < y_local.len() && global_j < y.len() {
443 y[global_j] += y_local[local_j];
444 }
445 }
446 }
447 }
448
449 pub fn get_cluster(&self, level: usize, index: usize) -> Option<&Cluster> {
451 self.cluster_levels.get(level)?.clusters.get(index)
452 }
453
454 pub fn num_clusters_at_level(&self, level: usize) -> usize {
456 self.cluster_levels
457 .get(level)
458 .map(|l| l.clusters.len())
459 .unwrap_or(0)
460 }
461}
462
463pub fn build_mlfmm_system(
471 elements: &[Element],
472 nodes: &Array2<f64>,
473 cluster_levels: Vec<ClusterLevel>,
474 physics: &PhysicsParams,
475) -> MlfmmSystem {
476 let num_dofs = count_dofs(elements);
477 let num_levels = cluster_levels.len();
478
479 let mut system = MlfmmSystem::new(num_dofs, num_levels, cluster_levels.clone());
480
481 if let Some(leaf_level) = cluster_levels.last() {
483 build_leaf_dof_mappings(&mut system, elements, &leaf_level.clusters);
484 }
485
486 if let Some(leaf_level) = cluster_levels.last() {
488 build_near_field(&mut system, elements, nodes, &leaf_level.clusters, physics);
489 }
490
491 for (level_idx, level) in cluster_levels.iter().enumerate().rev() {
493 let (sphere_coords, sphere_weights) =
494 unit_sphere_quadrature(level.theta_points, level.phi_points);
495
496 let t_matrices = build_t_matrices_level(
497 elements,
498 &level.clusters,
499 physics,
500 &sphere_coords,
501 &sphere_weights,
502 level_idx,
503 &cluster_levels,
504 );
505
506 system.t_matrices.push(t_matrices);
507 }
508 system.t_matrices.reverse(); for (level_idx, level) in cluster_levels.iter().enumerate() {
512 let (sphere_coords, _sphere_weights) =
513 unit_sphere_quadrature(level.theta_points, level.phi_points);
514
515 let d_matrices = build_d_matrices_level(
516 &level.clusters,
517 physics,
518 &sphere_coords,
519 level.expansion_terms,
520 level_idx,
521 );
522
523 system.d_matrices.push(d_matrices);
524 }
525
526 for (level_idx, level) in cluster_levels.iter().enumerate() {
528 let (sphere_coords, sphere_weights) =
529 unit_sphere_quadrature(level.theta_points, level.phi_points);
530
531 let s_matrices = build_s_matrices_level(
532 elements,
533 &level.clusters,
534 physics,
535 &sphere_coords,
536 &sphere_weights,
537 level_idx,
538 &cluster_levels,
539 );
540
541 system.s_matrices.push(s_matrices);
542 }
543
544 system
545}
546
547fn build_leaf_dof_mappings(
549 system: &mut MlfmmSystem,
550 elements: &[Element],
551 leaf_clusters: &[Cluster],
552) {
553 for cluster in leaf_clusters {
554 let mut dof_indices = Vec::new();
555 for &elem_idx in &cluster.element_indices {
556 let elem = &elements[elem_idx];
557 if elem.property.is_evaluation() {
558 continue;
559 }
560 dof_indices.extend(elem.dof_addresses.iter().copied());
562 }
563 system.leaf_dof_indices.push(dof_indices);
564 }
565}
566
567fn count_dofs(elements: &[Element]) -> usize {
569 elements
570 .iter()
571 .filter(|e| !e.property.is_evaluation())
572 .map(|e| e.dof_addresses.len())
573 .sum()
574}
575
576fn build_near_field(
578 system: &mut MlfmmSystem,
579 elements: &[Element],
580 nodes: &Array2<f64>,
581 leaf_clusters: &[Cluster],
582 physics: &PhysicsParams,
583) {
584 let gamma = Complex64::new(physics.gamma(), 0.0);
585 let tau = Complex64::new(physics.tau, 0.0);
586 let beta = physics.burton_miller_beta();
587
588 for (i, cluster_i) in leaf_clusters.iter().enumerate() {
590 let block = compute_near_block(
592 elements,
593 nodes,
594 &cluster_i.element_indices,
595 &cluster_i.element_indices,
596 physics,
597 gamma,
598 tau,
599 beta,
600 true,
601 );
602 system.near_matrix.push(NearFieldBlock {
603 source_cluster: i,
604 field_cluster: i,
605 coefficients: block,
606 });
607
608 for &j in &cluster_i.near_clusters {
610 if j > i {
611 let cluster_j = &leaf_clusters[j];
612 let block = compute_near_block(
613 elements,
614 nodes,
615 &cluster_i.element_indices,
616 &cluster_j.element_indices,
617 physics,
618 gamma,
619 tau,
620 beta,
621 false,
622 );
623 system.near_matrix.push(NearFieldBlock {
624 source_cluster: i,
625 field_cluster: j,
626 coefficients: block,
627 });
628 }
629 }
630 }
631}
632
633fn compute_near_block(
635 elements: &[Element],
636 nodes: &Array2<f64>,
637 source_indices: &[usize],
638 field_indices: &[usize],
639 physics: &PhysicsParams,
640 gamma: Complex64,
641 tau: Complex64,
642 beta: Complex64,
643 is_self: bool,
644) -> Array2<Complex64> {
645 let n_source = source_indices.len();
646 let n_field = field_indices.len();
647 let mut block = Array2::zeros((n_source, n_field));
648
649 for (i, &src_idx) in source_indices.iter().enumerate() {
650 let source_elem = &elements[src_idx];
651 if source_elem.property.is_evaluation() {
652 continue;
653 }
654
655 for (j, &fld_idx) in field_indices.iter().enumerate() {
656 let field_elem = &elements[fld_idx];
657 if field_elem.property.is_evaluation() {
658 continue;
659 }
660
661 let element_coords = get_element_coords(field_elem, nodes);
662
663 let result = if is_self && src_idx == fld_idx {
664 singular_integration(
666 &source_elem.center,
667 &source_elem.normal,
668 &element_coords,
669 field_elem.element_type,
670 physics,
671 None,
672 0,
673 false,
674 )
675 } else {
676 regular_integration(
678 &source_elem.center,
679 &source_elem.normal,
680 &element_coords,
681 field_elem.element_type,
682 field_elem.area,
683 physics,
684 None,
685 0,
686 false,
687 )
688 };
689
690 let coeff =
692 result.dg_dn_integral * gamma * tau + result.d2g_dnxdny_integral * beta;
693 block[[i, j]] = coeff;
694 }
695 }
696
697 block
698}
699
700fn build_t_matrices_level(
702 elements: &[Element],
703 clusters: &[Cluster],
704 physics: &PhysicsParams,
705 sphere_coords: &[[f64; 3]],
706 sphere_weights: &[f64],
707 level_idx: usize,
708 cluster_levels: &[ClusterLevel],
709) -> LevelMatrices {
710 let k = physics.wave_number;
711 let num_sphere_points = sphere_coords.len();
712 let num_clusters = clusters.len();
713 let mut matrices = Vec::with_capacity(num_clusters);
714
715 let is_leaf_level = level_idx == cluster_levels.len() - 1;
716
717 for cluster in clusters {
718 if is_leaf_level {
719 let num_elem = cluster.element_indices.len();
721 let mut t_matrix = Array2::zeros((num_sphere_points, num_elem));
722
723 for (j, &elem_idx) in cluster.element_indices.iter().enumerate() {
724 let elem = &elements[elem_idx];
725 if elem.property.is_evaluation() {
726 continue;
727 }
728
729 for (p, coord) in sphere_coords.iter().enumerate() {
730 let diff: Vec<f64> = (0..3)
731 .map(|d| elem.center[d] - cluster.center[d])
732 .collect();
733 let s_dot_diff: f64 = (0..3).map(|d| coord[d] * diff[d]).sum();
734
735 let exp_factor =
736 Complex64::new((k * s_dot_diff).cos(), -(k * s_dot_diff).sin());
737
738 t_matrix[[p, j]] = exp_factor * sphere_weights[p];
739 }
740 }
741
742 matrices.push(t_matrix);
743 } else {
744 let num_children = cluster.sons.len();
746 let child_num_points = if level_idx + 1 < cluster_levels.len() {
747 cluster_levels[level_idx + 1].theta_points
748 * cluster_levels[level_idx + 1].phi_points
749 } else {
750 num_sphere_points
751 };
752
753 let mut t_matrix = Array2::zeros((num_sphere_points, num_children * child_num_points));
754
755 for (child_idx, &child_cluster_idx) in cluster.sons.iter().enumerate() {
757 if level_idx + 1 < cluster_levels.len() {
758 let child_cluster = &cluster_levels[level_idx + 1].clusters[child_cluster_idx];
759
760 for (p, coord) in sphere_coords.iter().enumerate() {
762 let diff: Vec<f64> = (0..3)
763 .map(|d| child_cluster.center[d] - cluster.center[d])
764 .collect();
765 let s_dot_diff: f64 = (0..3).map(|d| coord[d] * diff[d]).sum();
766
767 let exp_factor =
768 Complex64::new((k * s_dot_diff).cos(), -(k * s_dot_diff).sin());
769
770 for cp in 0..child_num_points {
772 t_matrix[[p, child_idx * child_num_points + cp]] =
773 exp_factor * sphere_weights[p] / child_num_points as f64;
774 }
775 }
776 }
777 }
778
779 matrices.push(t_matrix);
780 }
781 }
782
783 LevelMatrices {
784 level: level_idx,
785 num_clusters,
786 matrices,
787 }
788}
789
790fn build_d_matrices_level(
792 clusters: &[Cluster],
793 physics: &PhysicsParams,
794 sphere_coords: &[[f64; 3]],
795 n_terms: usize,
796 level_idx: usize,
797) -> Vec<DMatrixEntry> {
798 let k = physics.wave_number;
799 let num_sphere_points = sphere_coords.len();
800 let mut d_entries = Vec::new();
801
802 for (i, cluster_i) in clusters.iter().enumerate() {
803 for &j in &cluster_i.far_clusters {
804 let cluster_j = &clusters[j];
805
806 let diff: Vec<f64> = (0..3)
808 .map(|d| cluster_i.center[d] - cluster_j.center[d])
809 .collect();
810 let r = (diff[0] * diff[0] + diff[1] * diff[1] + diff[2] * diff[2]).sqrt();
811
812 if r < 1e-15 {
813 continue; }
815
816 let kr = k * r;
817
818 let mut d_matrix = Array2::zeros((num_sphere_points, num_sphere_points));
819
820 let h_funcs = spherical_hankel_first_kind(n_terms.max(2), kr, 1.0);
822 let _p_funcs = legendre_polynomials(n_terms.max(2), 0.0);
823
824 for p in 0..num_sphere_points {
826 d_matrix[[p, p]] = h_funcs[0] * Complex64::new(0.0, k);
827 }
828
829 d_entries.push(DMatrixEntry {
830 source_cluster: i,
831 field_cluster: j,
832 level: level_idx,
833 coefficients: d_matrix,
834 });
835 }
836 }
837
838 d_entries
839}
840
841fn build_s_matrices_level(
843 elements: &[Element],
844 clusters: &[Cluster],
845 physics: &PhysicsParams,
846 sphere_coords: &[[f64; 3]],
847 sphere_weights: &[f64],
848 level_idx: usize,
849 cluster_levels: &[ClusterLevel],
850) -> LevelMatrices {
851 let k = physics.wave_number;
852 let num_sphere_points = sphere_coords.len();
853 let num_clusters = clusters.len();
854 let mut matrices = Vec::with_capacity(num_clusters);
855
856 let is_leaf_level = level_idx == cluster_levels.len() - 1;
857
858 for cluster in clusters {
859 if is_leaf_level {
860 let num_elem = cluster.element_indices.len();
862 let mut s_matrix = Array2::zeros((num_elem, num_sphere_points));
863
864 for (j, &elem_idx) in cluster.element_indices.iter().enumerate() {
865 let elem = &elements[elem_idx];
866 if elem.property.is_evaluation() {
867 continue;
868 }
869
870 for (p, coord) in sphere_coords.iter().enumerate() {
871 let diff: Vec<f64> = (0..3)
872 .map(|d| elem.center[d] - cluster.center[d])
873 .collect();
874 let s_dot_diff: f64 = (0..3).map(|d| coord[d] * diff[d]).sum();
875
876 let exp_factor =
877 Complex64::new((k * s_dot_diff).cos(), (k * s_dot_diff).sin());
878
879 s_matrix[[j, p]] = exp_factor * sphere_weights[p];
880 }
881 }
882
883 matrices.push(s_matrix);
884 } else {
885 let num_children = cluster.sons.len();
887 let child_num_points = if level_idx + 1 < cluster_levels.len() {
888 cluster_levels[level_idx + 1].theta_points
889 * cluster_levels[level_idx + 1].phi_points
890 } else {
891 num_sphere_points
892 };
893
894 let mut s_matrix = Array2::zeros((num_children * child_num_points, num_sphere_points));
895
896 for (child_idx, &child_cluster_idx) in cluster.sons.iter().enumerate() {
897 if level_idx + 1 < cluster_levels.len() {
898 let child_cluster = &cluster_levels[level_idx + 1].clusters[child_cluster_idx];
899
900 for (p, coord) in sphere_coords.iter().enumerate() {
901 let diff: Vec<f64> = (0..3)
902 .map(|d| child_cluster.center[d] - cluster.center[d])
903 .collect();
904 let s_dot_diff: f64 = (0..3).map(|d| coord[d] * diff[d]).sum();
905
906 let exp_factor =
907 Complex64::new((k * s_dot_diff).cos(), (k * s_dot_diff).sin());
908
909 for cp in 0..child_num_points {
910 s_matrix[[child_idx * child_num_points + cp, p]] =
911 exp_factor * sphere_weights[p] / child_num_points as f64;
912 }
913 }
914 }
915 }
916
917 matrices.push(s_matrix);
918 }
919 }
920
921 LevelMatrices {
922 level: level_idx,
923 num_clusters,
924 matrices,
925 }
926}
927
928fn get_element_coords(element: &Element, nodes: &Array2<f64>) -> Array2<f64> {
930 let num_nodes = element.connectivity.len();
931 let mut coords = Array2::zeros((num_nodes, 3));
932
933 for (i, &node_idx) in element.connectivity.iter().enumerate() {
934 for j in 0..3 {
935 coords[[i, j]] = nodes[[node_idx, j]];
936 }
937 }
938
939 coords
940}
941
942pub fn estimate_num_levels(
946 num_elements: usize,
947 elements_per_leaf: usize,
948 min_levels: usize,
949 max_levels: usize,
950) -> usize {
951 if num_elements == 0 {
952 return min_levels;
953 }
954
955 let mut levels = 1;
957 let mut n = num_elements;
958
959 while n > elements_per_leaf && levels < max_levels {
960 n /= 8;
961 levels += 1;
962 }
963
964 levels.max(min_levels).min(max_levels)
965}
966
967pub fn build_cluster_tree(
971 elements: &[Element],
972 target_elements_per_leaf: usize,
973 physics: &PhysicsParams,
974) -> Vec<ClusterLevel> {
975 let num_elements = elements.len();
976 let num_levels = estimate_num_levels(num_elements, target_elements_per_leaf, 1, 8);
977
978 let mut levels = Vec::with_capacity(num_levels);
979
980 let (min_coords, max_coords) = compute_bounding_box(elements);
982 let root_center = [
983 (min_coords[0] + max_coords[0]) / 2.0,
984 (min_coords[1] + max_coords[1]) / 2.0,
985 (min_coords[2] + max_coords[2]) / 2.0,
986 ];
987 let root_radius = ((max_coords[0] - min_coords[0]).powi(2)
988 + (max_coords[1] - min_coords[1]).powi(2)
989 + (max_coords[2] - min_coords[2]).powi(2))
990 .sqrt()
991 / 2.0;
992
993 let root_cluster = Cluster::new(Array1::from_vec(root_center.to_vec()));
995 let all_indices: Vec<usize> = (0..num_elements).collect();
996
997 let mut root_level = ClusterLevel::new(1);
999 let mut root = root_cluster;
1000 root.element_indices = all_indices;
1001 root.radius = root_radius;
1002 root.level = 0;
1003
1004 let kr = physics.wave_number * root_radius;
1006 root_level.expansion_terms = ((kr + 6.0 * kr.ln().max(1.0)) as usize).max(4).min(30);
1007 root_level.theta_points = root_level.expansion_terms;
1008 root_level.phi_points = 2 * root_level.expansion_terms;
1009
1010 root_level.clusters.push(root);
1011 root_level.num_original = 1;
1012 root_level.max_radius = root_radius;
1013 root_level.avg_radius = root_radius;
1014 root_level.min_radius = root_radius;
1015
1016 levels.push(root_level);
1017
1018 if num_levels > 1 {
1020 subdivide_level(&mut levels, elements, 0, target_elements_per_leaf, physics);
1021 }
1022
1023 for level in &mut levels {
1025 compute_near_far_lists(&mut level.clusters, physics.wave_number);
1026 }
1027
1028 levels
1029}
1030
1031fn compute_bounding_box(elements: &[Element]) -> ([f64; 3], [f64; 3]) {
1033 let mut min_coords = [f64::MAX, f64::MAX, f64::MAX];
1034 let mut max_coords = [f64::MIN, f64::MIN, f64::MIN];
1035
1036 for elem in elements {
1037 for d in 0..3 {
1038 min_coords[d] = min_coords[d].min(elem.center[d]);
1039 max_coords[d] = max_coords[d].max(elem.center[d]);
1040 }
1041 }
1042
1043 (min_coords, max_coords)
1044}
1045
1046fn subdivide_level(
1048 levels: &mut Vec<ClusterLevel>,
1049 elements: &[Element],
1050 parent_level: usize,
1051 target_elements_per_leaf: usize,
1052 physics: &PhysicsParams,
1053) {
1054 let parent_clusters = levels[parent_level].clusters.clone();
1055 let mut child_level = ClusterLevel::new(parent_clusters.len() * 8);
1056
1057 let mut max_radius = 0.0_f64;
1058 let mut min_radius = f64::MAX;
1059 let mut sum_radius = 0.0_f64;
1060
1061 for (parent_idx, parent) in parent_clusters.iter().enumerate() {
1062 if parent.element_indices.len() <= target_elements_per_leaf {
1063 let mut leaf = parent.clone();
1065 leaf.level = parent_level + 1;
1066 leaf.father = Some(parent_idx);
1067
1068 max_radius = max_radius.max(leaf.radius);
1069 min_radius = min_radius.min(leaf.radius);
1070 sum_radius += leaf.radius;
1071
1072 let child_idx = child_level.clusters.len();
1073 levels[parent_level].clusters[parent_idx].sons.push(child_idx);
1074 child_level.clusters.push(leaf);
1075 continue;
1076 }
1077
1078 let half_size = parent.radius / 2.0;
1080 let offsets = [
1081 [-1.0, -1.0, -1.0],
1082 [-1.0, -1.0, 1.0],
1083 [-1.0, 1.0, -1.0],
1084 [-1.0, 1.0, 1.0],
1085 [1.0, -1.0, -1.0],
1086 [1.0, -1.0, 1.0],
1087 [1.0, 1.0, -1.0],
1088 [1.0, 1.0, 1.0],
1089 ];
1090
1091 for offset in &offsets {
1092 let child_center = [
1093 parent.center[0] + offset[0] * half_size * 0.5,
1094 parent.center[1] + offset[1] * half_size * 0.5,
1095 parent.center[2] + offset[2] * half_size * 0.5,
1096 ];
1097
1098 let child_elements: Vec<usize> = parent
1100 .element_indices
1101 .iter()
1102 .filter(|&&idx| {
1103 let elem = &elements[idx];
1104 let dx = elem.center[0] - parent.center[0];
1105 let dy = elem.center[1] - parent.center[1];
1106 let dz = elem.center[2] - parent.center[2];
1107
1108 (dx * offset[0] >= 0.0)
1109 && (dy * offset[1] >= 0.0)
1110 && (dz * offset[2] >= 0.0)
1111 })
1112 .copied()
1113 .collect();
1114
1115 if child_elements.is_empty() {
1116 continue;
1117 }
1118
1119 let mut child = Cluster::new(Array1::from_vec(child_center.to_vec()));
1120 child.element_indices = child_elements;
1121 child.radius = half_size;
1122 child.level = parent_level + 1;
1123 child.father = Some(parent_idx);
1124
1125 max_radius = max_radius.max(child.radius);
1126 min_radius = min_radius.min(child.radius);
1127 sum_radius += child.radius;
1128
1129 let child_idx = child_level.clusters.len();
1130 levels[parent_level].clusters[parent_idx].sons.push(child_idx);
1131 child_level.clusters.push(child);
1132 }
1133 }
1134
1135 let num_clusters = child_level.clusters.len();
1136 child_level.num_original = num_clusters;
1137 child_level.max_radius = max_radius;
1138 child_level.min_radius = min_radius;
1139 child_level.avg_radius = if num_clusters > 0 {
1140 sum_radius / num_clusters as f64
1141 } else {
1142 0.0
1143 };
1144
1145 let kr = physics.wave_number * child_level.avg_radius;
1147 child_level.expansion_terms = ((kr + 6.0 * kr.ln().max(1.0)) as usize).max(4).min(30);
1148 child_level.theta_points = child_level.expansion_terms;
1149 child_level.phi_points = 2 * child_level.expansion_terms;
1150
1151 levels.push(child_level);
1152
1153 let current_level = levels.len() - 1;
1155 let should_continue = levels[current_level]
1156 .clusters
1157 .iter()
1158 .any(|c| c.element_indices.len() > target_elements_per_leaf);
1159
1160 if should_continue && current_level < 7 {
1161 subdivide_level(levels, elements, current_level, target_elements_per_leaf, physics);
1162 }
1163}
1164
1165fn compute_near_far_lists(clusters: &mut [Cluster], wave_number: f64) {
1167 let n = clusters.len();
1168 let min_separation = 2.0; for i in 0..n {
1171 let center_i = clusters[i].center.clone();
1172 let radius_i = clusters[i].radius;
1173
1174 let mut near = Vec::new();
1175 let mut far = Vec::new();
1176
1177 for j in 0..n {
1178 if i == j {
1179 continue;
1180 }
1181
1182 let center_j = &clusters[j].center;
1183 let radius_j = clusters[j].radius;
1184
1185 let dist = ((center_i[0] - center_j[0]).powi(2)
1186 + (center_i[1] - center_j[1]).powi(2)
1187 + (center_i[2] - center_j[2]).powi(2))
1188 .sqrt();
1189
1190 let separation = dist / (radius_i + radius_j).max(1e-15);
1191
1192 let kr = wave_number * dist;
1194 let is_acoustic_far = kr > 2.0;
1195
1196 if separation > min_separation && is_acoustic_far {
1197 far.push(j);
1198 } else {
1199 near.push(j);
1200 }
1201 }
1202
1203 clusters[i].near_clusters = near;
1204 clusters[i].far_clusters = far;
1205 }
1206}
1207
1208#[cfg(test)]
1209mod tests {
1210 use super::*;
1211 use crate::core::types::{BoundaryCondition, ElementProperty, ElementType};
1212 use ndarray::array;
1213
1214 fn make_test_elements() -> (Vec<Element>, Array2<f64>) {
1215 let nodes = Array2::from_shape_vec(
1216 (4, 3),
1217 vec![0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.5, 1.0, 0.0, 1.5, 1.0, 0.0],
1218 )
1219 .unwrap();
1220
1221 let elem0 = Element {
1222 connectivity: vec![0, 1, 2],
1223 element_type: ElementType::Tri3,
1224 property: ElementProperty::Surface,
1225 normal: array![0.0, 0.0, 1.0],
1226 node_normals: Array2::zeros((3, 3)),
1227 center: array![0.5, 1.0 / 3.0, 0.0],
1228 area: 0.5,
1229 boundary_condition: BoundaryCondition::Velocity(vec![Complex64::new(1.0, 0.0)]),
1230 group: 0,
1231 dof_addresses: vec![0],
1232 };
1233
1234 let elem1 = Element {
1235 connectivity: vec![1, 3, 2],
1236 element_type: ElementType::Tri3,
1237 property: ElementProperty::Surface,
1238 normal: array![0.0, 0.0, 1.0],
1239 node_normals: Array2::zeros((3, 3)),
1240 center: array![1.0, 2.0 / 3.0, 0.0],
1241 area: 0.5,
1242 boundary_condition: BoundaryCondition::Velocity(vec![Complex64::new(0.0, 0.0)]),
1243 group: 0,
1244 dof_addresses: vec![1],
1245 };
1246
1247 (vec![elem0, elem1], nodes)
1248 }
1249
1250 #[test]
1251 fn test_estimate_num_levels() {
1252 assert_eq!(estimate_num_levels(10, 10, 1, 8), 1);
1253 assert_eq!(estimate_num_levels(100, 10, 1, 8), 3);
1255 assert!(estimate_num_levels(1000, 10, 1, 8) >= 3);
1257 }
1258
1259 #[test]
1260 fn test_build_cluster_tree() {
1261 let (elements, _nodes) = make_test_elements();
1262 let physics = PhysicsParams::new(100.0, 343.0, 1.21, false);
1263
1264 let tree = build_cluster_tree(&elements, 10, &physics);
1265
1266 assert!(!tree.is_empty());
1267 assert!(!tree[0].clusters.is_empty());
1268 assert_eq!(tree[0].clusters[0].element_indices.len(), elements.len());
1269 }
1270
1271 #[test]
1272 fn test_build_mlfmm_system() {
1273 let (elements, nodes) = make_test_elements();
1274 let physics = PhysicsParams::new(100.0, 343.0, 1.21, false);
1275
1276 let cluster_tree = build_cluster_tree(&elements, 10, &physics);
1277 let system = build_mlfmm_system(&elements, &nodes, cluster_tree, &physics);
1278
1279 assert_eq!(system.num_dofs, 2);
1280 assert!(system.num_levels >= 1);
1281 }
1282
1283 #[test]
1284 fn test_mlfmm_matvec() {
1285 let (elements, nodes) = make_test_elements();
1286 let physics = PhysicsParams::new(100.0, 343.0, 1.21, false);
1287
1288 let cluster_tree = build_cluster_tree(&elements, 10, &physics);
1289 let system = build_mlfmm_system(&elements, &nodes, cluster_tree, &physics);
1290
1291 let x = Array1::from_vec(vec![Complex64::new(1.0, 0.0), Complex64::new(0.0, 1.0)]);
1292 let y = system.matvec(&x);
1293
1294 assert_eq!(y.len(), 2);
1295 }
1296}