bem/core/postprocess/
pressure.rs

1//! Field pressure evaluation at exterior points
2//!
3//! Computes the acoustic pressure field at evaluation points using
4//! the BEM representation formula:
5//!
6//! p(x) = p_inc(x) + ∫_Γ [p(y) ∂G/∂n_y - (∂p/∂n)(y) G(x,y)] dS_y
7//!
8//! For a rigid scatterer (Neumann BC), the normal velocity is zero,
9//! so the formula simplifies to:
10//!
11//! p(x) = p_inc(x) + ∫_Γ p(y) ∂G/∂n_y dS_y
12
13use ndarray::{Array1, Array2};
14use num_complex::Complex64;
15use std::f64::consts::PI;
16
17// Note: Using inline Green's function computation for performance
18use crate::core::incident::IncidentField;
19use crate::core::integration::gauss::triangle_quadrature;
20use crate::core::types::{Element, ElementType, PhysicsParams};
21
22/// Field evaluation result at a single point
23#[derive(Debug, Clone)]
24pub struct FieldPoint {
25    /// Position coordinates [x, y, z]
26    pub position: [f64; 3],
27    /// Incident pressure
28    pub p_incident: Complex64,
29    /// Scattered pressure
30    pub p_scattered: Complex64,
31    /// Total pressure (incident + scattered)
32    pub p_total: Complex64,
33}
34
35impl FieldPoint {
36    /// Create new field point
37    pub fn new(position: [f64; 3], p_incident: Complex64, p_scattered: Complex64) -> Self {
38        Self {
39            position,
40            p_incident,
41            p_scattered,
42            p_total: p_incident + p_scattered,
43        }
44    }
45
46    /// Get pressure magnitude in dB SPL (re: 20 μPa)
47    pub fn spl_db(&self) -> f64 {
48        let p_ref = 20e-6; // Reference pressure for dB SPL
49        20.0 * (self.p_total.norm() / p_ref).log10()
50    }
51
52    /// Get pressure magnitude (absolute value)
53    pub fn magnitude(&self) -> f64 {
54        self.p_total.norm()
55    }
56
57    /// Get pressure phase in radians
58    pub fn phase(&self) -> f64 {
59        self.p_total.arg()
60    }
61}
62
63/// Compute scattered field at evaluation points using surface solution
64///
65/// Uses the Kirchhoff-Helmholtz integral equation for exterior problems:
66/// p_scat(x) = ∫_Γ [p(y) ∂G/∂n_y(x,y) - v_n(y) G(x,y)] dS_y
67///
68/// This implementation uses proper Gauss quadrature over each element,
69/// matching the C++ NumCalc NC_IntegrationTBEM function for accuracy.
70///
71/// # Arguments
72/// * `eval_points` - Evaluation points (N × 3 array)
73/// * `elements` - Mesh elements with surface solution
74/// * `nodes` - Node coordinates
75/// * `surface_pressure` - Solved surface pressure values (one per element)
76/// * `surface_velocity` - Surface normal velocity (usually from BC, can be zero for rigid)
77/// * `physics` - Physical parameters
78///
79/// # Returns
80/// Scattered pressure at each evaluation point
81pub fn compute_scattered_field(
82    eval_points: &Array2<f64>,
83    elements: &[Element],
84    nodes: &Array2<f64>,
85    surface_pressure: &Array1<Complex64>,
86    surface_velocity: Option<&Array1<Complex64>>,
87    physics: &PhysicsParams,
88) -> Array1<Complex64> {
89    let n_eval = eval_points.nrows();
90    let k = physics.wave_number;
91    let harmonic_factor = physics.harmonic_factor;
92    let wavruim = k * harmonic_factor;
93
94    let mut p_scattered = Array1::zeros(n_eval);
95
96    // Get boundary elements (exclude evaluation-only elements)
97    let boundary_elements: Vec<_> = elements
98        .iter()
99        .filter(|e| !e.property.is_evaluation())
100        .collect();
101
102    for i in 0..n_eval {
103        let x = Array1::from_vec(vec![
104            eval_points[[i, 0]],
105            eval_points[[i, 1]],
106            eval_points[[i, 2]],
107        ]);
108
109        let mut p_scat = Complex64::new(0.0, 0.0);
110
111        for (j, elem) in boundary_elements.iter().enumerate() {
112            // Get surface values (constant over element)
113            let p_surf = surface_pressure[j];
114            let v_surf = surface_velocity.map_or(Complex64::new(0.0, 0.0), |v| v[j]);
115
116            // Get element node coordinates
117            let elem_coords = get_element_coords(elem, nodes);
118
119            // Integrate using Gauss quadrature (matching NumCalc approach)
120            let contribution = integrate_element_field(
121                &x,
122                &elem_coords,
123                elem.element_type,
124                p_surf,
125                v_surf,
126                k,
127                wavruim,
128            );
129
130            p_scat += contribution;
131        }
132
133        p_scattered[i] = p_scat;
134    }
135
136    p_scattered
137}
138
139/// Get element node coordinates from connectivity
140fn get_element_coords(elem: &Element, nodes: &Array2<f64>) -> Array2<f64> {
141    let n_nodes = elem.connectivity.len();
142    let mut coords = Array2::zeros((n_nodes, 3));
143    for (i, &node_idx) in elem.connectivity.iter().enumerate() {
144        for j in 0..3 {
145            coords[[i, j]] = nodes[[node_idx, j]];
146        }
147    }
148    coords
149}
150
151/// Integrate element contribution using Gauss quadrature
152///
153/// This matches the NC_IntegrationTBEM function from NumCalc.
154fn integrate_element_field(
155    x: &Array1<f64>,
156    elem_coords: &Array2<f64>,
157    elem_type: ElementType,
158    p_surf: Complex64,
159    v_surf: Complex64,
160    _k: f64,
161    wavruim: f64,
162) -> Complex64 {
163    let mut result = Complex64::new(0.0, 0.0);
164
165    // Use order 3 quadrature (7-point for triangles) for good accuracy
166    let gauss_order = 3;
167    let quad_points = match elem_type {
168        ElementType::Tri3 => triangle_quadrature(gauss_order),
169        ElementType::Quad4 => {
170            // For quads, use tensor product rule (approximate)
171            triangle_quadrature(gauss_order)
172        }
173    };
174
175    for (xi, eta, weight) in quad_points {
176        // Compute shape functions and derivatives
177        let (shape_fn, shape_ds, shape_dt) = match elem_type {
178            ElementType::Tri3 => {
179                let shape_fn = vec![1.0 - xi - eta, xi, eta];
180                let shape_ds = vec![-1.0, 1.0, 0.0];
181                let shape_dt = vec![-1.0, 0.0, 1.0];
182                (shape_fn, shape_ds, shape_dt)
183            }
184            ElementType::Quad4 => {
185                // Use triangle approximation for quads
186                let shape_fn = vec![1.0 - xi - eta, xi, eta];
187                let shape_ds = vec![-1.0, 1.0, 0.0];
188                let shape_dt = vec![-1.0, 0.0, 1.0];
189                (shape_fn, shape_ds, shape_dt)
190            }
191        };
192
193        // Compute position at quadrature point
194        let mut crd_poi: Array1<f64> = Array1::zeros(3);
195        let mut dx_ds: Array1<f64> = Array1::zeros(3);
196        let mut dx_dt: Array1<f64> = Array1::zeros(3);
197
198        let n_nodes = elem_coords.nrows().min(3); // Limit to 3 for tri3
199        for n in 0..n_nodes {
200            for d in 0..3 {
201                crd_poi[d] += shape_fn[n] * elem_coords[[n, d]];
202                dx_ds[d] += shape_ds[n] * elem_coords[[n, d]];
203                dx_dt[d] += shape_dt[n] * elem_coords[[n, d]];
204            }
205        }
206
207        // Compute normal and Jacobian
208        let normal: Array1<f64> = Array1::from_vec(vec![
209            dx_ds[1] * dx_dt[2] - dx_ds[2] * dx_dt[1],
210            dx_ds[2] * dx_dt[0] - dx_ds[0] * dx_dt[2],
211            dx_ds[0] * dx_dt[1] - dx_ds[1] * dx_dt[0],
212        ]);
213        let jacobian = normal.dot(&normal).sqrt();
214
215        if jacobian < 1e-15 {
216            continue;
217        }
218
219        let el_norm = &normal / jacobian;
220
221        // Distance from evaluation point to quadrature point
222        let mut r_vec: Array1<f64> = Array1::zeros(3);
223        for d in 0..3 {
224            r_vec[d] = crd_poi[d] - x[d];
225        }
226        let r = r_vec.dot(&r_vec).sqrt();
227
228        if r < 1e-15 {
229            continue;
230        }
231
232        // Weight factor
233        let vjacwe = jacobian * weight;
234
235        // Green's function: G = exp(ikr) / (4πr)
236        let kr = wavruim * r;
237        let re1 = 4.0 * PI * r;
238        let zgrfu = Complex64::new(kr.cos() / re1, kr.sin() / re1);
239
240        // Derivative factor: z1 = -1/r + ik
241        let z1 = Complex64::new(-1.0 / r, wavruim);
242        let zgikr = zgrfu * z1;
243
244        // ∂r/∂n_y = (y-x)·n_y / r
245        let drdn_y = r_vec.dot(&el_norm) / r;
246
247        // Double layer contribution: p * ∂G/∂n_y * dS
248        let zdgrdn = zgikr * drdn_y;
249        result += p_surf * zdgrdn * vjacwe;
250
251        // Single layer contribution: -v * G * dS
252        // (v_surf is normally zero for rigid scatterer)
253        if v_surf.norm() > 1e-15 {
254            result -= v_surf * zgrfu * vjacwe;
255        }
256    }
257
258    result
259}
260
261/// Compute total field (incident + scattered) at evaluation points
262///
263/// # Arguments
264/// * `eval_points` - Evaluation points (N × 3 array)
265/// * `elements` - Mesh elements
266/// * `nodes` - Node coordinates
267/// * `surface_pressure` - Solved surface pressure
268/// * `incident_field` - Incident field specification
269/// * `physics` - Physical parameters
270///
271/// # Returns
272/// Vector of FieldPoint with all pressure components
273pub fn compute_total_field(
274    eval_points: &Array2<f64>,
275    elements: &[Element],
276    nodes: &Array2<f64>,
277    surface_pressure: &Array1<Complex64>,
278    surface_velocity: Option<&Array1<Complex64>>,
279    incident_field: &IncidentField,
280    physics: &PhysicsParams,
281) -> Vec<FieldPoint> {
282    // Compute incident field
283    let p_incident = incident_field.evaluate_pressure(eval_points, physics);
284
285    // Compute scattered field
286    let p_scattered = compute_scattered_field(
287        eval_points,
288        elements,
289        nodes,
290        surface_pressure,
291        surface_velocity,
292        physics,
293    );
294
295    // Combine results
296    let n_eval = eval_points.nrows();
297    let mut results = Vec::with_capacity(n_eval);
298
299    for i in 0..n_eval {
300        let position = [
301            eval_points[[i, 0]],
302            eval_points[[i, 1]],
303            eval_points[[i, 2]],
304        ];
305        results.push(FieldPoint::new(position, p_incident[i], p_scattered[i]));
306    }
307
308    results
309}
310
311/// Generate evaluation points on a sphere around the origin
312///
313/// # Arguments
314/// * `radius` - Sphere radius
315/// * `n_theta` - Number of polar divisions
316/// * `n_phi` - Number of azimuthal divisions
317///
318/// # Returns
319/// Array of points on the sphere (N × 3)
320pub fn generate_sphere_eval_points(radius: f64, n_theta: usize, n_phi: usize) -> Array2<f64> {
321    let mut points = Vec::new();
322
323    for i in 0..n_theta {
324        let theta = PI * (i as f64 + 0.5) / n_theta as f64;
325        let sin_theta = theta.sin();
326        let cos_theta = theta.cos();
327
328        for j in 0..n_phi {
329            let phi = 2.0 * PI * j as f64 / n_phi as f64;
330
331            points.push(radius * sin_theta * phi.cos());
332            points.push(radius * sin_theta * phi.sin());
333            points.push(radius * cos_theta);
334        }
335    }
336
337    let n_points = n_theta * n_phi;
338    Array2::from_shape_vec((n_points, 3), points).unwrap()
339}
340
341/// Generate evaluation points along a line
342///
343/// # Arguments
344/// * `start` - Start point
345/// * `end` - End point
346/// * `n_points` - Number of points
347///
348/// # Returns
349/// Array of points along the line (N × 3)
350pub fn generate_line_eval_points(start: [f64; 3], end: [f64; 3], n_points: usize) -> Array2<f64> {
351    let mut points = Vec::new();
352
353    for i in 0..n_points {
354        let t = i as f64 / (n_points - 1).max(1) as f64;
355        points.push(start[0] + t * (end[0] - start[0]));
356        points.push(start[1] + t * (end[1] - start[1]));
357        points.push(start[2] + t * (end[2] - start[2]));
358    }
359
360    Array2::from_shape_vec((n_points, 3), points).unwrap()
361}
362
363/// Generate evaluation points on a plane (for field maps)
364///
365/// # Arguments
366/// * `center` - Center of the plane
367/// * `normal` - Normal to the plane (determines orientation)
368/// * `extent` - Half-size of the plane in each direction
369/// * `n_points` - Number of points in each direction
370///
371/// # Returns
372/// Array of points on the plane (N² × 3)
373pub fn generate_plane_eval_points(
374    center: [f64; 3],
375    normal: [f64; 3],
376    extent: f64,
377    n_points: usize,
378) -> Array2<f64> {
379    // Find two vectors perpendicular to normal
380    let n = Array1::from_vec(vec![normal[0], normal[1], normal[2]]);
381    let n_norm = n.dot(&n).sqrt();
382    let n = &n / n_norm;
383
384    // Choose an arbitrary vector not parallel to n
385    let arbitrary = if n[0].abs() < 0.9 {
386        Array1::from_vec(vec![1.0, 0.0, 0.0])
387    } else {
388        Array1::from_vec(vec![0.0, 1.0, 0.0])
389    };
390
391    // First basis vector (perpendicular to n)
392    let u = {
393        let cross = Array1::from_vec(vec![
394            n[1] * arbitrary[2] - n[2] * arbitrary[1],
395            n[2] * arbitrary[0] - n[0] * arbitrary[2],
396            n[0] * arbitrary[1] - n[1] * arbitrary[0],
397        ]);
398        let norm = cross.dot(&cross).sqrt();
399        cross / norm
400    };
401
402    // Second basis vector (perpendicular to both n and u)
403    let v = Array1::from_vec(vec![
404        n[1] * u[2] - n[2] * u[1],
405        n[2] * u[0] - n[0] * u[2],
406        n[0] * u[1] - n[1] * u[0],
407    ]);
408
409    let mut points = Vec::new();
410
411    for i in 0..n_points {
412        let s = -extent + 2.0 * extent * i as f64 / (n_points - 1).max(1) as f64;
413        for j in 0..n_points {
414            let t = -extent + 2.0 * extent * j as f64 / (n_points - 1).max(1) as f64;
415
416            points.push(center[0] + s * u[0] + t * v[0]);
417            points.push(center[1] + s * u[1] + t * v[1]);
418            points.push(center[2] + s * u[2] + t * v[2]);
419        }
420    }
421
422    Array2::from_shape_vec((n_points * n_points, 3), points).unwrap()
423}
424
425/// Compute RCS (Radar Cross Section) from far-field pattern
426///
427/// For acoustic scattering, the RCS is defined as:
428/// σ = lim_{r→∞} 4πr² |p_scat/p_inc|²
429///
430/// # Arguments
431/// * `surface_pressure` - Solved surface pressure
432/// * `elements` - Mesh elements
433/// * `direction` - Far-field direction (unit vector)
434/// * `physics` - Physical parameters
435///
436/// # Returns
437/// RCS value
438pub fn compute_rcs(
439    surface_pressure: &Array1<Complex64>,
440    elements: &[Element],
441    direction: [f64; 3],
442    physics: &PhysicsParams,
443) -> f64 {
444    let k = physics.wave_number;
445
446    // Far-field pattern is computed using the stationary phase approximation
447    // For constant elements: F(θ,φ) = Σ p_j * exp(-ik r_j · d) * A_j * (n_j · d)
448
449    let mut far_field = Complex64::new(0.0, 0.0);
450
451    let boundary_elements: Vec<_> = elements
452        .iter()
453        .filter(|e| !e.property.is_evaluation())
454        .collect();
455
456    let d = Array1::from_vec(vec![direction[0], direction[1], direction[2]]);
457
458    for (j, elem) in boundary_elements.iter().enumerate() {
459        let center = &elem.center;
460        let normal = &elem.normal;
461        let area = elem.area;
462
463        // Phase term: exp(-ik r · d)
464        let phase = -k * center.dot(&d);
465        let exp_phase = Complex64::new(phase.cos(), phase.sin());
466
467        // Normal factor: n · d (double layer contribution)
468        let n_dot_d = normal.dot(&d);
469
470        // Contribution
471        let ik = Complex64::new(0.0, k);
472        far_field += surface_pressure[j] * exp_phase * area * ik * n_dot_d;
473    }
474
475    // RCS = |F|² / |p_inc|²
476    // For unit amplitude incident wave, RCS = |F|²
477    4.0 * PI * far_field.norm_sqr()
478}
479
480#[cfg(test)]
481mod tests {
482    use super::*;
483    use crate::core::types::ElementProperty;
484    use ndarray::array;
485
486    fn make_physics(k: f64) -> PhysicsParams {
487        let c = 343.0;
488        let freq = k * c / (2.0 * PI);
489        PhysicsParams::new(freq, c, 1.21, false)
490    }
491
492    #[test]
493    fn test_field_point_creation() {
494        let p_inc = Complex64::new(1.0, 0.0);
495        let p_scat = Complex64::new(0.5, 0.3);
496
497        let fp = FieldPoint::new([1.0, 0.0, 0.0], p_inc, p_scat);
498
499        assert!((fp.p_total - (p_inc + p_scat)).norm() < 1e-10);
500        assert!(fp.magnitude() > 0.0);
501    }
502
503    #[test]
504    fn test_generate_sphere_eval_points() {
505        let points = generate_sphere_eval_points(2.0, 10, 20);
506
507        assert_eq!(points.nrows(), 200);
508        assert_eq!(points.ncols(), 3);
509
510        // Check all points are at radius 2.0
511        for i in 0..points.nrows() {
512            let r =
513                (points[[i, 0]].powi(2) + points[[i, 1]].powi(2) + points[[i, 2]].powi(2)).sqrt();
514            assert!((r - 2.0).abs() < 1e-10);
515        }
516    }
517
518    #[test]
519    fn test_generate_line_eval_points() {
520        let points = generate_line_eval_points([0.0, 0.0, 0.0], [1.0, 0.0, 0.0], 11);
521
522        assert_eq!(points.nrows(), 11);
523
524        // First point should be at origin
525        assert!(points[[0, 0]].abs() < 1e-10);
526
527        // Last point should be at (1, 0, 0)
528        assert!((points[[10, 0]] - 1.0).abs() < 1e-10);
529
530        // Point spacing should be uniform
531        let spacing = points[[1, 0]] - points[[0, 0]];
532        assert!((spacing - 0.1).abs() < 1e-10);
533    }
534
535    #[test]
536    fn test_generate_plane_eval_points() {
537        let points = generate_plane_eval_points([0.0, 0.0, 0.0], [0.0, 0.0, 1.0], 1.0, 5);
538
539        assert_eq!(points.nrows(), 25);
540        assert_eq!(points.ncols(), 3);
541
542        // All points should have z = 0 (on the xy plane)
543        for i in 0..points.nrows() {
544            assert!(points[[i, 2]].abs() < 1e-10);
545        }
546    }
547
548    #[test]
549    fn test_scattered_field_basic() {
550        // Simple test with a single element
551        use crate::core::types::BoundaryCondition;
552
553        let physics = make_physics(1.0);
554
555        // Single triangular element
556        let nodes =
557            Array2::from_shape_vec((3, 3), vec![0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.5, 1.0, 0.0])
558                .unwrap();
559
560        let elem = Element {
561            connectivity: vec![0, 1, 2],
562            element_type: crate::core::types::ElementType::Tri3,
563            property: ElementProperty::Surface,
564            normal: array![0.0, 0.0, 1.0],
565            node_normals: Array2::zeros((3, 3)),
566            center: array![0.5, 1.0 / 3.0, 0.0],
567            area: 0.5,
568            boundary_condition: BoundaryCondition::Velocity(vec![Complex64::new(0.0, 0.0)]),
569            group: 0,
570            dof_addresses: vec![0],
571        };
572
573        let elements = vec![elem];
574        let surface_pressure = Array1::from_vec(vec![Complex64::new(1.0, 0.0)]);
575
576        // Evaluation point above the element
577        let eval_points = Array2::from_shape_vec((1, 3), vec![0.5, 0.5, 1.0]).unwrap();
578
579        let p_scat = compute_scattered_field(
580            &eval_points,
581            &elements,
582            &nodes,
583            &surface_pressure,
584            None,
585            &physics,
586        );
587
588        // Scattered field should be non-zero
589        assert!(p_scat[0].norm() > 0.0);
590    }
591}