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, ElementType, 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    // Loop over source elements
109    for (iel, source_elem) in elements.iter().enumerate() {
110        // Skip evaluation elements (property == 2)
111        if source_elem.property.is_evaluation() {
112            continue;
113        }
114
115        let source_point = &source_elem.center;
116        let source_normal = &source_elem.normal;
117        let source_dof = source_elem.dof_addresses[0];
118
119        // Get boundary condition type and value
120        let (bc_type, bc_value) = get_bc_type_and_value(&source_elem.boundary_condition);
121
122        // Add diagonal free terms (jump terms from BIE)
123        add_free_terms(
124            &mut system,
125            source_dof,
126            bc_type,
127            &bc_value,
128            gamma,
129            tau,
130            beta,
131        );
132
133        // Loop over field elements
134        for (jel, field_elem) in elements.iter().enumerate() {
135            // Skip evaluation elements
136            if field_elem.property.is_evaluation() {
137                continue;
138            }
139
140            // Get field element coordinates
141            let element_coords = get_element_coords(field_elem, nodes);
142            let field_dof = field_elem.dof_addresses[0];
143
144            // Get field BC for RHS contribution
145            let (field_bc_type, field_bc_values) = get_bc_type_and_value(&field_elem.boundary_condition);
146            let compute_rhs = has_nonzero_bc(&field_bc_values);
147
148            // Compute integrals
149            let result = if jel == iel {
150                // Singular integration (self-element)
151                singular_integration(
152                    source_point,
153                    source_normal,
154                    &element_coords,
155                    field_elem.element_type,
156                    physics,
157                    if compute_rhs { Some(&field_bc_values) } else { None },
158                    field_bc_type,
159                    compute_rhs,
160                )
161            } else {
162                // Regular integration
163                regular_integration(
164                    source_point,
165                    source_normal,
166                    &element_coords,
167                    field_elem.element_type,
168                    field_elem.area,
169                    physics,
170                    if compute_rhs { Some(&field_bc_values) } else { None },
171                    field_bc_type,
172                    compute_rhs,
173                )
174            };
175
176            // Assemble contributions to matrix and RHS
177            assemble_tbem(
178                &mut system,
179                source_dof,
180                field_dof,
181                bc_type,
182                field_bc_type,
183                &result,
184                gamma,
185                tau,
186                beta,
187                compute_rhs,
188            );
189        }
190    }
191
192    system
193}
194
195/// Count total number of DOFs
196fn count_dofs(elements: &[Element]) -> usize {
197    elements
198        .iter()
199        .filter(|e| !e.property.is_evaluation())
200        .map(|e| e.dof_addresses.len())
201        .sum()
202}
203
204/// Get boundary condition type (0=velocity, 1=pressure) and values
205fn get_bc_type_and_value(bc: &BoundaryCondition) -> (i32, Vec<Complex64>) {
206    match bc {
207        BoundaryCondition::Velocity(v) => (0, v.clone()),
208        BoundaryCondition::Pressure(p) => (1, p.clone()),
209        BoundaryCondition::VelocityWithAdmittance { velocity, .. } => (0, velocity.clone()),
210        BoundaryCondition::TransferAdmittance { .. } => (2, vec![Complex64::new(0.0, 0.0)]),
211        BoundaryCondition::TransferWithSurfaceAdmittance { .. } => (2, vec![Complex64::new(0.0, 0.0)]),
212    }
213}
214
215/// Check if BC values are non-zero
216fn has_nonzero_bc(values: &[Complex64]) -> bool {
217    values.iter().any(|v| v.norm() > 1e-15)
218}
219
220/// Get element node coordinates as Array2
221fn get_element_coords(element: &Element, nodes: &Array2<f64>) -> Array2<f64> {
222    let num_nodes = element.connectivity.len();
223    let mut coords = Array2::zeros((num_nodes, 3));
224
225    for (i, &node_idx) in element.connectivity.iter().enumerate() {
226        for j in 0..3 {
227            coords[[i, j]] = nodes[[node_idx, j]];
228        }
229    }
230
231    coords
232}
233
234/// Add diagonal free terms from jump conditions
235///
236/// From NC_ComputeEntriesForTBEM in C++:
237/// - For velocity BC: diagonal += (Admia3*zBta3 - Gama3)*0.5
238///   Without admittance: diagonal += -γ/2
239/// - For pressure BC: diagonal -= β*τ/2
240///
241/// This matches the classic Burton-Miller formulation for exterior problems.
242fn add_free_terms(
243    system: &mut TbemSystem,
244    source_dof: usize,
245    bc_type: i32,
246    bc_value: &[Complex64],
247    gamma: Complex64,
248    tau: Complex64,
249    beta: Complex64,
250) {
251    let avg_bc = bc_value.iter().sum::<Complex64>() / bc_value.len() as f64;
252
253    match bc_type {
254        0 => {
255            // Velocity BC: diagonal term from CBIE jump
256            // Using c = -0.5 (matching C++ NumCalc) with negated K' gives best results
257            system.matrix[[source_dof, source_dof]] -= gamma * 0.5;
258            // RHS contribution from velocity BC (from NC_ComputeEntriesForTBEM line 239)
259            // z0 += Zbvi03*zBta3*Tao_*0.5
260            system.rhs[source_dof] += avg_bc * beta * tau * 0.5;
261        }
262        1 => {
263            // Pressure BC: diagonal term from HBIE jump
264            // C++: zcoefl[diagonal] -= zBta3*(0.5*Tao_)
265            system.matrix[[source_dof, source_dof]] -= beta * tau * 0.5;
266            // RHS contribution from pressure BC
267            system.rhs[source_dof] += avg_bc * tau * 0.5;
268        }
269        _ => {
270            // Transfer admittance - more complex handling
271        }
272    }
273}
274
275/// Assemble element contributions to matrix and RHS
276///
277/// The Burton-Miller formulation gives:
278/// - For velocity BC on field element: coeff += (γτH + βE)
279/// - For pressure BC on field element: coeff += (-γτG - βH^T)
280fn assemble_tbem(
281    system: &mut TbemSystem,
282    source_dof: usize,
283    field_dof: usize,
284    _source_bc_type: i32,
285    field_bc_type: i32,
286    result: &IntegrationResult,
287    gamma: Complex64,
288    tau: Complex64,
289    beta: Complex64,
290    compute_rhs: bool,
291) {
292    let coeff = match field_bc_type {
293        0 => {
294            // Velocity BC (Neumann): unknown is total surface pressure p
295            // DIRECT formulation using Kirchhoff-Helmholtz representation
296            // Burton-Miller: (c + K' + βH)[p] = p_inc + β*∂p_inc/∂n
297            // K' = ∫ ∂G/∂n_y dS = dg_dn_integral (double-layer transpose)
298            // H = ∫ ∂²G/(∂n_x∂n_y) dS = d2g_dnxdny_integral (hypersingular)
299            -result.dg_dn_integral * gamma * tau + result.d2g_dnxdny_integral * beta
300        }
301        1 => {
302            // Pressure BC (Dirichlet): unknown is surface velocity v
303            // Uses G and H^T kernels
304            -(result.g_integral * gamma * tau + result.dg_dnx_integral * beta)
305        }
306        _ => Complex64::new(0.0, 0.0),
307    };
308
309    system.matrix[[source_dof, field_dof]] += coeff;
310
311    if compute_rhs {
312        system.rhs[source_dof] += result.rhs_contribution;
313    }
314}
315
316/// Build TBEM system in parallel using rayon
317#[cfg(feature = "parallel")]
318pub fn build_tbem_system_parallel(
319    elements: &[Element],
320    nodes: &Array2<f64>,
321    physics: &PhysicsParams,
322) -> TbemSystem {
323    use rayon::prelude::*;
324    use std::sync::Mutex;
325
326    let num_dofs = count_dofs(elements);
327    let system = Mutex::new(TbemSystem::new(num_dofs));
328
329    let gamma = Complex64::new(physics.gamma(), 0.0);
330    let tau = Complex64::new(physics.tau, 0.0);
331    let beta = physics.burton_miller_beta();
332
333    // Process source elements in parallel
334    elements.par_iter().enumerate().for_each(|(iel, source_elem)| {
335        if source_elem.property.is_evaluation() {
336            return;
337        }
338
339        let source_point = &source_elem.center;
340        let source_normal = &source_elem.normal;
341        let source_dof = source_elem.dof_addresses[0];
342        let (bc_type, bc_value) = get_bc_type_and_value(&source_elem.boundary_condition);
343
344        // Compute row locally
345        let mut local_row = Array1::<Complex64>::zeros(num_dofs);
346        let mut local_rhs = Complex64::new(0.0, 0.0);
347
348        // Add free terms (c = +1/2 for exterior problem)
349        let avg_bc = bc_value.iter().sum::<Complex64>() / bc_value.len() as f64;
350        match bc_type {
351            0 => {
352                local_row[source_dof] += gamma * 0.5;
353                local_rhs += avg_bc * beta * tau * 0.5;
354            }
355            1 => {
356                local_row[source_dof] += beta * tau * 0.5;
357                local_rhs += avg_bc * tau * 0.5;
358            }
359            _ => {}
360        }
361
362        // Loop over field elements
363        for (jel, field_elem) in elements.iter().enumerate() {
364            if field_elem.property.is_evaluation() {
365                continue;
366            }
367
368            let element_coords = get_element_coords(field_elem, nodes);
369            let field_dof = field_elem.dof_addresses[0];
370            let (field_bc_type, field_bc_values) = get_bc_type_and_value(&field_elem.boundary_condition);
371            let compute_rhs = has_nonzero_bc(&field_bc_values);
372
373            let result = if jel == iel {
374                singular_integration(
375                    source_point,
376                    source_normal,
377                    &element_coords,
378                    field_elem.element_type,
379                    physics,
380                    if compute_rhs { Some(&field_bc_values) } else { None },
381                    field_bc_type,
382                    compute_rhs,
383                )
384            } else {
385                regular_integration(
386                    source_point,
387                    source_normal,
388                    &element_coords,
389                    field_elem.element_type,
390                    field_elem.area,
391                    physics,
392                    if compute_rhs { Some(&field_bc_values) } else { None },
393                    field_bc_type,
394                    compute_rhs,
395                )
396            };
397
398            let coeff = match field_bc_type {
399                // Velocity BC: K' + βH for direct formulation
400                0 => result.dg_dn_integral * gamma * tau + result.d2g_dnxdny_integral * beta,
401                // Pressure BC: G + βH^T for direct formulation
402                1 => -(result.g_integral * gamma * tau + result.dg_dnx_integral * beta),
403                _ => Complex64::new(0.0, 0.0),
404            };
405
406            local_row[field_dof] += coeff;
407
408            if compute_rhs {
409                local_rhs += result.rhs_contribution;
410            }
411        }
412
413        // Write row to global system
414        let mut sys = system.lock().unwrap();
415        for j in 0..num_dofs {
416            sys.matrix[[source_dof, j]] += local_row[j];
417        }
418        sys.rhs[source_dof] += local_rhs;
419    });
420
421    system.into_inner().unwrap()
422}
423
424/// Apply row sum correction to enforce the theoretical row sum property
425///
426/// For a closed surface with exterior Burton-Miller formulation:
427/// - Row sum should be approximately 0 (c + K[1] + β*E[1] ≈ 0.5 - 0.5 + 0 = 0)
428///
429/// Due to numerical integration errors in E[1], the computed row sum may be nonzero.
430/// This correction subtracts the row sum from the diagonal to enforce the property.
431///
432/// # Arguments
433/// * `system` - The assembled TBEM system (modified in place)
434///
435/// # Returns
436/// The average row sum magnitude before correction (for diagnostics)
437pub fn apply_row_sum_correction(system: &mut TbemSystem) -> f64 {
438    let n = system.num_dofs;
439    let mut max_row_sum_norm = 0.0f64;
440    let mut total_row_sum = Complex64::new(0.0, 0.0);
441
442    for i in 0..n {
443        // Compute row sum
444        let mut row_sum = Complex64::new(0.0, 0.0);
445        for j in 0..n {
446            row_sum += system.matrix[[i, j]];
447        }
448
449        total_row_sum += row_sum;
450        max_row_sum_norm = max_row_sum_norm.max(row_sum.norm());
451
452        // Correct diagonal to make row sum zero
453        system.matrix[[i, i]] -= row_sum;
454    }
455
456    let avg_row_sum = total_row_sum.norm() / n as f64;
457    avg_row_sum
458}
459
460/// Build TBEM system with row sum correction
461///
462/// This is recommended for closed surface scattering problems where
463/// the E[1] = 0 property should hold but numerical errors cause drift.
464pub fn build_tbem_system_corrected(
465    elements: &[Element],
466    nodes: &Array2<f64>,
467    physics: &PhysicsParams,
468) -> (TbemSystem, f64) {
469    let mut system = build_tbem_system(elements, nodes, physics);
470    let correction = apply_row_sum_correction(&mut system);
471    (system, correction)
472}
473
474#[cfg(test)]
475mod tests {
476    use super::*;
477    use crate::core::types::ElementProperty;
478    use ndarray::array;
479
480    fn make_simple_mesh() -> (Vec<Element>, Array2<f64>) {
481        // Simple 2-element mesh for testing
482        let nodes = Array2::from_shape_vec(
483            (4, 3),
484            vec![
485                0.0, 0.0, 0.0,  // node 0
486                1.0, 0.0, 0.0,  // node 1
487                0.5, 1.0, 0.0,  // node 2
488                1.5, 1.0, 0.0,  // node 3
489            ],
490        )
491        .unwrap();
492
493        let elem0 = Element {
494            connectivity: vec![0, 1, 2],
495            element_type: ElementType::Tri3,
496            property: ElementProperty::Surface,
497            normal: array![0.0, 0.0, 1.0],
498            node_normals: Array2::zeros((3, 3)),
499            center: array![0.5, 1.0 / 3.0, 0.0],
500            area: 0.5,
501            boundary_condition: BoundaryCondition::Velocity(vec![Complex64::new(1.0, 0.0)]),
502            group: 0,
503            dof_addresses: vec![0],
504        };
505
506        let elem1 = Element {
507            connectivity: vec![1, 3, 2],
508            element_type: ElementType::Tri3,
509            property: ElementProperty::Surface,
510            normal: array![0.0, 0.0, 1.0],
511            node_normals: Array2::zeros((3, 3)),
512            center: array![1.0, 2.0 / 3.0, 0.0],
513            area: 0.5,
514            boundary_condition: BoundaryCondition::Velocity(vec![Complex64::new(0.0, 0.0)]),
515            group: 0,
516            dof_addresses: vec![1],
517        };
518
519        (vec![elem0, elem1], nodes)
520    }
521
522    #[test]
523    fn test_build_tbem_system() {
524        let (elements, nodes) = make_simple_mesh();
525        let physics = PhysicsParams::new(100.0, 343.0, 1.21, false);
526
527        let system = build_tbem_system(&elements, &nodes, &physics);
528
529        assert_eq!(system.num_dofs, 2);
530        assert_eq!(system.matrix.shape(), &[2, 2]);
531        assert_eq!(system.rhs.len(), 2);
532
533        // Matrix should have non-zero diagonal entries
534        assert!(system.matrix[[0, 0]].norm() > 1e-15);
535        assert!(system.matrix[[1, 1]].norm() > 1e-15);
536    }
537
538    #[test]
539    fn test_count_dofs() {
540        let (elements, _) = make_simple_mesh();
541        assert_eq!(count_dofs(&elements), 2);
542    }
543
544    #[test]
545    fn test_get_element_coords() {
546        let (elements, nodes) = make_simple_mesh();
547        let coords = get_element_coords(&elements[0], &nodes);
548
549        assert_eq!(coords.shape(), &[3, 3]);
550        assert!((coords[[0, 0]] - 0.0).abs() < 1e-10); // node 0
551        assert!((coords[[1, 0]] - 1.0).abs() < 1e-10); // node 1
552    }
553}