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