Skip to main content

math_audio_bem/core/assembly/
tbem.rs

1//! Traditional BEM (O(N²) dense matrix) assembly
2//!
3//! Direct port of NC_BuildSystemTBEM from NC_EquationSystem.cpp.
4//! Uses the Burton-Miller formulation to avoid spurious resonances.
5
6use ndarray::{Array1, Array2};
7use num_complex::Complex64;
8
9use crate::core::integration::{regular_integration, singular_integration};
10use crate::core::types::{BoundaryCondition, Element, IntegrationResult, PhysicsParams};
11
12/// Result of TBEM assembly
13pub struct TbemSystem {
14    /// Coefficient matrix (dense, row-major)
15    pub matrix: Array2<Complex64>,
16    /// Right-hand side vector
17    pub rhs: Array1<Complex64>,
18    /// Number of DOFs
19    pub num_dofs: usize,
20}
21
22impl TbemSystem {
23    /// Create a new system with given number of DOFs
24    pub fn new(num_dofs: usize) -> Self {
25        Self {
26            matrix: Array2::zeros((num_dofs, num_dofs)),
27            rhs: Array1::zeros(num_dofs),
28            num_dofs,
29        }
30    }
31}
32
33/// Build the TBEM system matrix and RHS vector
34///
35/// This implements NC_BuildSystemTBEM from the C++ code.
36/// Uses Burton-Miller formulation: (γτG + βH^T) for the system matrix.
37///
38/// # Arguments
39/// * `elements` - Vector of mesh elements
40/// * `nodes` - Node coordinates (num_nodes × 3)
41/// * `physics` - Physics parameters (wave number, etc.)
42///
43/// # Returns
44/// TbemSystem with assembled matrix and RHS
45pub fn build_tbem_system(
46    elements: &[Element],
47    nodes: &Array2<f64>,
48    physics: &PhysicsParams,
49) -> TbemSystem {
50    build_tbem_system_with_beta(elements, nodes, physics, physics.burton_miller_beta())
51}
52
53/// Build TBEM system with mesh-adaptive Burton-Miller coupling
54///
55/// Uses β = i/(k + k_ref) where k_ref = 1/element_size for better conditioning
56/// at low-to-mid frequencies (ka < 2). This avoids the 1/k divergence of the
57/// traditional β = i/k formulation.
58///
59/// # Arguments
60/// * `elements` - Vector of mesh elements
61/// * `nodes` - Node coordinates (num_nodes × 3)
62/// * `physics` - Physics parameters (wave number, etc.)
63/// * `avg_element_size` - Average element edge length or sqrt(area)
64pub fn build_tbem_system_bounded(
65    elements: &[Element],
66    nodes: &Array2<f64>,
67    physics: &PhysicsParams,
68    avg_element_size: f64,
69) -> TbemSystem {
70    let beta = physics.burton_miller_beta_optimal(avg_element_size);
71    build_tbem_system_with_beta(elements, nodes, physics, beta)
72}
73
74/// Build TBEM system with scaled Burton-Miller coupling
75///
76/// Uses β = scale × i/k for improved numerical conditioning.
77/// Testing shows scale = 2.0 gives significantly better accuracy than the
78/// standard β = i/k, especially at low to medium frequencies (ka < 2).
79///
80/// # Arguments
81/// * `elements` - Vector of mesh elements
82/// * `nodes` - Node coordinates (num_nodes × 3)
83/// * `physics` - Physics parameters (wave number, etc.)
84/// * `scale` - Beta scaling factor (recommended: 2.0)
85pub fn build_tbem_system_scaled(
86    elements: &[Element],
87    nodes: &Array2<f64>,
88    physics: &PhysicsParams,
89    scale: f64,
90) -> TbemSystem {
91    let beta = physics.burton_miller_beta_scaled(scale);
92    build_tbem_system_with_beta(elements, nodes, physics, beta)
93}
94
95/// Build TBEM system with custom Burton-Miller coupling parameter
96pub fn build_tbem_system_with_beta(
97    elements: &[Element],
98    nodes: &Array2<f64>,
99    physics: &PhysicsParams,
100    beta: Complex64,
101) -> TbemSystem {
102    let num_dofs = count_dofs(elements);
103    let mut system = TbemSystem::new(num_dofs);
104
105    let gamma = Complex64::new(physics.gamma(), 0.0);
106    let tau = Complex64::new(physics.tau, 0.0);
107
108    // Calculate characteristic radius from mesh bounding box
109    // This is more robust than averaging element center distances
110    let mut min_coord = [f64::MAX, f64::MAX, f64::MAX];
111    let mut max_coord = [f64::MIN, f64::MIN, f64::MIN];
112
113    for element in elements.iter() {
114        for d in 0..3 {
115            min_coord[d] = min_coord[d].min(element.center[d]);
116            max_coord[d] = max_coord[d].max(element.center[d]);
117        }
118    }
119
120    // Use half the maximum extent as characteristic radius
121    let characteristic_radius = 0.5
122        * (max_coord[0] - min_coord[0])
123            .max(max_coord[1] - min_coord[1])
124            .max(max_coord[2] - min_coord[2]);
125
126    let ka = physics.wave_number * characteristic_radius;
127
128    // Switch between formulations based on frequency
129    // Low freq (ka < 0.5): Use Modified formulation (+K) for stability
130    // High freq (ka >= 0.5): Use Standard formulation (-K) for accuracy
131    let dg_dn_sign = if ka < 0.5 { 1.0 } else { -1.0 };
132
133    // Loop over source elements
134    for (iel, source_elem) in elements.iter().enumerate() {
135        // Skip evaluation elements (property == 2)
136        if source_elem.property.is_evaluation() {
137            continue;
138        }
139
140        let source_point = &source_elem.center;
141        let source_normal = &source_elem.normal;
142        let source_dof = source_elem.dof_addresses[0];
143
144        // Get boundary condition type and value
145        let (bc_type, bc_value) = get_bc_type_and_value(&source_elem.boundary_condition);
146
147        // Add diagonal free terms (jump terms from BIE)
148        add_free_terms(
149            &mut system,
150            source_dof,
151            bc_type,
152            &bc_value,
153            gamma,
154            tau,
155            beta,
156        );
157
158        // Loop over field elements
159        for (jel, field_elem) in elements.iter().enumerate() {
160            // Skip evaluation elements
161            if field_elem.property.is_evaluation() {
162                continue;
163            }
164
165            // Get field element coordinates
166            let element_coords = get_element_coords(field_elem, nodes);
167            let field_dof = field_elem.dof_addresses[0];
168
169            // Get field BC for RHS contribution
170            let (field_bc_type, field_bc_values) =
171                get_bc_type_and_value(&field_elem.boundary_condition);
172            let compute_rhs = has_nonzero_bc(&field_bc_values);
173
174            // Compute integrals
175            let mut result = if jel == iel {
176                // Singular integration (self-element)
177                singular_integration(
178                    source_point,
179                    source_normal,
180                    &element_coords,
181                    field_elem.element_type,
182                    physics,
183                    if compute_rhs {
184                        Some(&field_bc_values)
185                    } else {
186                        None
187                    },
188                    field_bc_type,
189                    compute_rhs,
190                )
191            } else {
192                // Regular integration
193                regular_integration(
194                    source_point,
195                    source_normal,
196                    &element_coords,
197                    field_elem.element_type,
198                    field_elem.area,
199                    physics,
200                    if compute_rhs {
201                        Some(&field_bc_values)
202                    } else {
203                        None
204                    },
205                    field_bc_type,
206                    compute_rhs,
207                )
208            };
209
210            // Apply sign switch for stability/accuracy trade-off
211            result.dg_dn_integral *= dg_dn_sign;
212
213            // Assemble contributions to matrix and RHS
214            assemble_tbem(
215                &mut system,
216                source_dof,
217                field_dof,
218                bc_type,
219                field_bc_type,
220                &result,
221                gamma,
222                tau,
223                beta,
224                compute_rhs,
225            );
226        }
227    }
228
229    system
230}
231
232/// Count total number of DOFs
233fn count_dofs(elements: &[Element]) -> usize {
234    elements
235        .iter()
236        .filter(|e| !e.property.is_evaluation())
237        .map(|e| e.dof_addresses.len())
238        .sum()
239}
240
241/// Get boundary condition type (0=velocity, 1=pressure) and values
242fn get_bc_type_and_value(bc: &BoundaryCondition) -> (i32, Vec<Complex64>) {
243    match bc {
244        BoundaryCondition::Velocity(v) => (0, v.clone()),
245        BoundaryCondition::Pressure(p) => (1, p.clone()),
246        BoundaryCondition::VelocityWithAdmittance { velocity, .. } => (0, velocity.clone()),
247        BoundaryCondition::TransferAdmittance { .. } => (2, vec![Complex64::new(0.0, 0.0)]),
248        BoundaryCondition::TransferWithSurfaceAdmittance { .. } => {
249            (2, vec![Complex64::new(0.0, 0.0)])
250        }
251    }
252}
253
254/// Check if BC values are non-zero
255fn has_nonzero_bc(values: &[Complex64]) -> bool {
256    values.iter().any(|v| v.norm() > 1e-15)
257}
258
259/// Get element node coordinates as Array2
260fn get_element_coords(element: &Element, nodes: &Array2<f64>) -> Array2<f64> {
261    let num_nodes = element.connectivity.len();
262    let mut coords = Array2::zeros((num_nodes, 3));
263
264    for (i, &node_idx) in element.connectivity.iter().enumerate() {
265        for j in 0..3 {
266            coords[[i, j]] = nodes[[node_idx, j]];
267        }
268    }
269
270    coords
271}
272
273/// Add diagonal free terms from jump conditions
274///
275/// From NC_ComputeEntriesForTBEM in C++:
276/// - For velocity BC: diagonal += (Admia3*zBta3 - Gama3)*0.5
277///   Without admittance: diagonal += -γ/2
278/// - For pressure BC: diagonal -= β*τ/2
279///
280/// This matches the classic Burton-Miller formulation for exterior problems.
281fn add_free_terms(
282    system: &mut TbemSystem,
283    source_dof: usize,
284    bc_type: i32,
285    bc_value: &[Complex64],
286    gamma: Complex64,
287    tau: Complex64,
288    beta: Complex64,
289) {
290    let avg_bc = bc_value.iter().sum::<Complex64>() / bc_value.len() as f64;
291
292    match bc_type {
293        0 => {
294            // Velocity BC: diagonal term from CBIE jump
295            // Using c = -0.5 (matching C++ NumCalc) with negated K' gives best results
296            system.matrix[[source_dof, source_dof]] -= gamma * 0.5;
297            // RHS contribution from velocity BC (from NC_ComputeEntriesForTBEM line 239)
298            // z0 += Zbvi03*zBta3*Tao_*0.5
299            system.rhs[source_dof] += avg_bc * beta * tau * 0.5;
300        }
301        1 => {
302            // Pressure BC: diagonal term from HBIE jump
303            // C++: zcoefl[diagonal] -= zBta3*(0.5*Tao_)
304            system.matrix[[source_dof, source_dof]] -= beta * tau * 0.5;
305            // RHS contribution from pressure BC
306            system.rhs[source_dof] += avg_bc * tau * 0.5;
307        }
308        _ => {
309            unimplemented!("Transfer admittance boundary conditions are not yet supported");
310        }
311    }
312}
313
314/// Assemble element contributions to matrix and RHS
315///
316/// The Burton-Miller formulation gives:
317/// - For velocity BC on field element: coeff += (γτH + βE)
318/// - For pressure BC on field element: coeff += (-γτG - βH^T)
319fn assemble_tbem(
320    system: &mut TbemSystem,
321    source_dof: usize,
322    field_dof: usize,
323    _source_bc_type: i32,
324    field_bc_type: i32,
325    result: &IntegrationResult,
326    gamma: Complex64,
327    tau: Complex64,
328    beta: Complex64,
329    compute_rhs: bool,
330) {
331    let coeff = match field_bc_type {
332        0 => {
333            // Velocity BC (Neumann): unknown is total surface pressure p
334            // DIRECT formulation using Kirchhoff-Helmholtz representation
335            // Burton-Miller: (c + K' + βH)[p] = p_inc + β*∂p_inc/∂n
336            // K' = ∫ ∂G/∂n_y dS = dg_dn_integral (double-layer transpose)
337            // H = ∫ ∂²G/(∂n_x∂n_y) dS = d2g_dnxdny_integral (hypersingular)
338            result.dg_dn_integral * gamma * tau + result.d2g_dnxdny_integral * beta
339        }
340        1 => {
341            // Pressure BC (Dirichlet): unknown is surface velocity v
342            // Uses G and H^T kernels
343            -(result.g_integral * gamma * tau + result.dg_dnx_integral * beta)
344        }
345        _ => Complex64::new(0.0, 0.0),
346    };
347
348    system.matrix[[source_dof, field_dof]] += coeff;
349
350    if compute_rhs {
351        system.rhs[source_dof] += result.rhs_contribution;
352    }
353}
354
355/// Build TBEM system in parallel using rayon
356#[cfg(feature = "parallel")]
357pub fn build_tbem_system_parallel(
358    elements: &[Element],
359    nodes: &Array2<f64>,
360    physics: &PhysicsParams,
361) -> TbemSystem {
362    use rayon::prelude::*;
363    use std::sync::Mutex;
364
365    let num_dofs = count_dofs(elements);
366    let system = Mutex::new(TbemSystem::new(num_dofs));
367
368    let gamma = Complex64::new(physics.gamma(), 0.0);
369    let tau = Complex64::new(physics.tau, 0.0);
370    let beta = physics.burton_miller_beta();
371
372    // Calculate characteristic radius from mesh bounding box (same as sequential)
373    let mut min_coord = [f64::MAX, f64::MAX, f64::MAX];
374    let mut max_coord = [f64::MIN, f64::MIN, f64::MIN];
375
376    for element in elements.iter() {
377        for d in 0..3 {
378            min_coord[d] = min_coord[d].min(element.center[d]);
379            max_coord[d] = max_coord[d].max(element.center[d]);
380        }
381    }
382
383    let characteristic_radius = 0.5
384        * (max_coord[0] - min_coord[0])
385            .max(max_coord[1] - min_coord[1])
386            .max(max_coord[2] - min_coord[2]);
387
388    let ka = physics.wave_number * characteristic_radius;
389
390    // Switch between formulations based on frequency
391    let dg_dn_sign = if ka < 0.5 { 1.0 } else { -1.0 };
392
393    // Process source elements in parallel
394    elements
395        .par_iter()
396        .enumerate()
397        .for_each(|(iel, source_elem)| {
398            if source_elem.property.is_evaluation() {
399                return;
400            }
401
402            let source_point = &source_elem.center;
403            let source_normal = &source_elem.normal;
404            let source_dof = source_elem.dof_addresses[0];
405            let (bc_type, bc_value) = get_bc_type_and_value(&source_elem.boundary_condition);
406
407            // Compute row locally
408            let mut local_row = Array1::<Complex64>::zeros(num_dofs);
409            let mut local_rhs = Complex64::new(0.0, 0.0);
410
411            // Add free terms (c = -0.5 for exterior problem with CBIE jump condition)
412            let avg_bc = bc_value.iter().sum::<Complex64>() / bc_value.len() as f64;
413            match bc_type {
414                0 => {
415                    local_row[source_dof] -= gamma * 0.5;
416                    local_rhs += avg_bc * beta * tau * 0.5;
417                }
418                1 => {
419                    local_row[source_dof] -= beta * tau * 0.5;
420                    local_rhs += avg_bc * tau * 0.5;
421                }
422                _ => {}
423            }
424
425            // Loop over field elements
426            for (jel, field_elem) in elements.iter().enumerate() {
427                if field_elem.property.is_evaluation() {
428                    continue;
429                }
430
431                let element_coords = get_element_coords(field_elem, nodes);
432                let field_dof = field_elem.dof_addresses[0];
433                let (field_bc_type, field_bc_values) =
434                    get_bc_type_and_value(&field_elem.boundary_condition);
435                let compute_rhs = has_nonzero_bc(&field_bc_values);
436
437                let mut result = if jel == iel {
438                    singular_integration(
439                        source_point,
440                        source_normal,
441                        &element_coords,
442                        field_elem.element_type,
443                        physics,
444                        if compute_rhs {
445                            Some(&field_bc_values)
446                        } else {
447                            None
448                        },
449                        field_bc_type,
450                        compute_rhs,
451                    )
452                } else {
453                    regular_integration(
454                        source_point,
455                        source_normal,
456                        &element_coords,
457                        field_elem.element_type,
458                        field_elem.area,
459                        physics,
460                        if compute_rhs {
461                            Some(&field_bc_values)
462                        } else {
463                            None
464                        },
465                        field_bc_type,
466                        compute_rhs,
467                    )
468                };
469
470                // Apply sign switch
471                result.dg_dn_integral *= dg_dn_sign;
472
473                let coeff = match field_bc_type {
474                    // Velocity BC: K' + βH for direct formulation
475                    0 => result.dg_dn_integral * gamma * tau + result.d2g_dnxdny_integral * beta,
476                    // Pressure BC: G + βH^T for direct formulation
477                    1 => -(result.g_integral * gamma * tau + result.dg_dnx_integral * beta),
478                    _ => Complex64::new(0.0, 0.0),
479                };
480
481                local_row[field_dof] += coeff;
482
483                if compute_rhs {
484                    local_rhs += result.rhs_contribution;
485                }
486            }
487
488            // Write row to global system
489            let mut sys = system.lock().unwrap();
490            for j in 0..num_dofs {
491                sys.matrix[[source_dof, j]] += local_row[j];
492            }
493            sys.rhs[source_dof] += local_rhs;
494        });
495
496    system.into_inner().unwrap()
497}
498
499/// Apply row sum correction to enforce the theoretical row sum property
500///
501/// For a closed surface with exterior Burton-Miller formulation:
502/// - Row sum should be approximately 0 (c + K\[1\] + β*E\[1\] ≈ 0.5 - 0.5 + 0 = 0)
503///
504/// Due to numerical integration errors in E\[1\], the computed row sum may be nonzero.
505/// This correction subtracts the row sum from the diagonal to enforce the property.
506///
507/// # Arguments
508/// * `system` - The assembled TBEM system (modified in place)
509///
510/// # Returns
511/// The average row sum magnitude before correction (for diagnostics)
512pub fn apply_row_sum_correction(system: &mut TbemSystem) -> f64 {
513    let n = system.num_dofs;
514    let mut max_row_sum_norm = 0.0f64;
515    let mut total_row_sum = Complex64::new(0.0, 0.0);
516
517    for i in 0..n {
518        // Compute row sum
519        let mut row_sum = Complex64::new(0.0, 0.0);
520        for j in 0..n {
521            row_sum += system.matrix[[i, j]];
522        }
523
524        total_row_sum += row_sum;
525        max_row_sum_norm = max_row_sum_norm.max(row_sum.norm());
526
527        // Correct diagonal to make row sum zero
528        system.matrix[[i, i]] -= row_sum;
529    }
530
531    total_row_sum.norm() / n as f64
532}
533
534/// Build TBEM system with row sum correction
535///
536/// This is recommended for closed surface scattering problems where
537/// the E\[1\] = 0 property should hold but numerical errors cause drift.
538pub fn build_tbem_system_corrected(
539    elements: &[Element],
540    nodes: &Array2<f64>,
541    physics: &PhysicsParams,
542) -> (TbemSystem, f64) {
543    let mut system = build_tbem_system(elements, nodes, physics);
544    let correction = apply_row_sum_correction(&mut system);
545    (system, correction)
546}
547
548#[cfg(test)]
549mod tests {
550    use super::*;
551    use crate::core::types::{ElementProperty, ElementType};
552    use ndarray::array;
553
554    fn make_simple_mesh() -> (Vec<Element>, Array2<f64>) {
555        // Simple 2-element mesh for testing
556        let nodes = Array2::from_shape_vec(
557            (4, 3),
558            vec![
559                0.0, 0.0, 0.0, // node 0
560                1.0, 0.0, 0.0, // node 1
561                0.5, 1.0, 0.0, // node 2
562                1.5, 1.0, 0.0, // node 3
563            ],
564        )
565        .unwrap();
566
567        let elem0 = Element {
568            connectivity: vec![0, 1, 2],
569            element_type: ElementType::Tri3,
570            property: ElementProperty::Surface,
571            normal: array![0.0, 0.0, 1.0],
572            node_normals: Array2::zeros((3, 3)),
573            center: array![0.5, 1.0 / 3.0, 0.0],
574            area: 0.5,
575            boundary_condition: BoundaryCondition::Velocity(vec![Complex64::new(1.0, 0.0)]),
576            group: 0,
577            dof_addresses: vec![0],
578        };
579
580        let elem1 = Element {
581            connectivity: vec![1, 3, 2],
582            element_type: ElementType::Tri3,
583            property: ElementProperty::Surface,
584            normal: array![0.0, 0.0, 1.0],
585            node_normals: Array2::zeros((3, 3)),
586            center: array![1.0, 2.0 / 3.0, 0.0],
587            area: 0.5,
588            boundary_condition: BoundaryCondition::Velocity(vec![Complex64::new(0.0, 0.0)]),
589            group: 0,
590            dof_addresses: vec![1],
591        };
592
593        (vec![elem0, elem1], nodes)
594    }
595
596    #[test]
597    fn test_build_tbem_system() {
598        let (elements, nodes) = make_simple_mesh();
599        let physics = PhysicsParams::new(100.0, 343.0, 1.21, false);
600
601        let system = build_tbem_system(&elements, &nodes, &physics);
602
603        assert_eq!(system.num_dofs, 2);
604        assert_eq!(system.matrix.shape(), &[2, 2]);
605        assert_eq!(system.rhs.len(), 2);
606
607        // Matrix should have non-zero diagonal entries
608        assert!(system.matrix[[0, 0]].norm() > 1e-15);
609        assert!(system.matrix[[1, 1]].norm() > 1e-15);
610    }
611
612    #[test]
613    fn test_count_dofs() {
614        let (elements, _) = make_simple_mesh();
615        assert_eq!(count_dofs(&elements), 2);
616    }
617
618    #[test]
619    fn test_get_element_coords() {
620        let (elements, nodes) = make_simple_mesh();
621        let coords = get_element_coords(&elements[0], &nodes);
622
623        assert_eq!(coords.shape(), &[3, 3]);
624        assert!((coords[[0, 0]] - 0.0).abs() < 1e-10); // node 0
625        assert!((coords[[1, 0]] - 1.0).abs() < 1e-10); // node 1
626    }
627}