bem/core/
incident.rs

1//! Incident field computation
2//!
3//! Computes incident acoustic fields for BEM excitation sources:
4//! - Plane waves
5//! - Point sources (monopoles)
6//!
7//! These are used to form the right-hand side of the BEM system.
8
9use ndarray::{Array1, Array2};
10use num_complex::Complex64;
11use std::f64::consts::PI;
12
13use crate::core::types::PhysicsParams;
14
15/// Incident field source type
16#[derive(Debug, Clone)]
17pub enum IncidentField {
18    /// Plane wave: p = A * exp(i k·x)
19    PlaneWave {
20        /// Direction of propagation (unit vector)
21        direction: [f64; 3],
22        /// Complex amplitude
23        amplitude: Complex64,
24    },
25
26    /// Point source (monopole): p = A * exp(ikr) / (4πr)
27    PointSource {
28        /// Source position
29        position: [f64; 3],
30        /// Source strength (volume velocity amplitude)
31        strength: Complex64,
32    },
33
34    /// Multiple plane waves
35    MultiplePlaneWaves(Vec<([f64; 3], Complex64)>),
36
37    /// Multiple point sources
38    MultiplePointSources(Vec<([f64; 3], Complex64)>),
39}
40
41impl IncidentField {
42    /// Create a plane wave with unit amplitude traveling in +z direction
43    ///
44    /// This matches the standard convention for Mie scattering theory where
45    /// the incident wave travels toward +z (positive z direction).
46    pub fn plane_wave_z() -> Self {
47        IncidentField::PlaneWave {
48            direction: [0.0, 0.0, 1.0],
49            amplitude: Complex64::new(1.0, 0.0),
50        }
51    }
52
53    /// Create a plane wave with unit amplitude traveling in -z direction
54    pub fn plane_wave_neg_z() -> Self {
55        IncidentField::PlaneWave {
56            direction: [0.0, 0.0, -1.0],
57            amplitude: Complex64::new(1.0, 0.0),
58        }
59    }
60
61    /// Create a plane wave with specified direction and amplitude
62    pub fn plane_wave(direction: [f64; 3], amplitude: f64) -> Self {
63        // Normalize direction
64        let len = (direction[0].powi(2) + direction[1].powi(2) + direction[2].powi(2)).sqrt();
65        let dir = if len > 1e-10 {
66            [direction[0] / len, direction[1] / len, direction[2] / len]
67        } else {
68            [0.0, 0.0, -1.0]
69        };
70
71        IncidentField::PlaneWave {
72            direction: dir,
73            amplitude: Complex64::new(amplitude, 0.0),
74        }
75    }
76
77    /// Create a point source at given position
78    pub fn point_source(position: [f64; 3], strength: f64) -> Self {
79        IncidentField::PointSource {
80            position,
81            strength: Complex64::new(strength, 0.0),
82        }
83    }
84
85    /// Evaluate incident pressure at given points
86    ///
87    /// # Arguments
88    /// * `points` - Evaluation points (N × 3 array)
89    /// * `physics` - Physical parameters (contains wave number k)
90    ///
91    /// # Returns
92    /// Complex pressure values at each point
93    pub fn evaluate_pressure(
94        &self,
95        points: &Array2<f64>,
96        physics: &PhysicsParams,
97    ) -> Array1<Complex64> {
98        let n = points.nrows();
99        let k = physics.wave_number;
100        let mut pressure = Array1::zeros(n);
101
102        match self {
103            IncidentField::PlaneWave { direction, amplitude } => {
104                for i in 0..n {
105                    // k·x = k * (d·x)
106                    let kdotx = k * (direction[0] * points[[i, 0]]
107                        + direction[1] * points[[i, 1]]
108                        + direction[2] * points[[i, 2]]);
109
110                    // p = A * exp(i k·x)
111                    pressure[i] = *amplitude * Complex64::new(kdotx.cos(), kdotx.sin());
112                }
113            }
114
115            IncidentField::PointSource { position, strength } => {
116                for i in 0..n {
117                    let dx = points[[i, 0]] - position[0];
118                    let dy = points[[i, 1]] - position[1];
119                    let dz = points[[i, 2]] - position[2];
120                    let r = (dx * dx + dy * dy + dz * dz).sqrt();
121
122                    if r > 1e-10 {
123                        let kr = k * r;
124                        // p = S * exp(ikr) / (4πr)
125                        let green = Complex64::new(kr.cos(), kr.sin()) / (4.0 * PI * r);
126                        pressure[i] = *strength * green;
127                    }
128                }
129            }
130
131            IncidentField::MultiplePlaneWaves(waves) => {
132                for (direction, amplitude) in waves {
133                    for i in 0..n {
134                        let kdotx = k * (direction[0] * points[[i, 0]]
135                            + direction[1] * points[[i, 1]]
136                            + direction[2] * points[[i, 2]]);
137                        pressure[i] += *amplitude * Complex64::new(kdotx.cos(), kdotx.sin());
138                    }
139                }
140            }
141
142            IncidentField::MultiplePointSources(sources) => {
143                for (position, strength) in sources {
144                    for i in 0..n {
145                        let dx = points[[i, 0]] - position[0];
146                        let dy = points[[i, 1]] - position[1];
147                        let dz = points[[i, 2]] - position[2];
148                        let r = (dx * dx + dy * dy + dz * dz).sqrt();
149
150                        if r > 1e-10 {
151                            let kr = k * r;
152                            let green = Complex64::new(kr.cos(), kr.sin()) / (4.0 * PI * r);
153                            pressure[i] += *strength * green;
154                        }
155                    }
156                }
157            }
158        }
159
160        pressure
161    }
162
163    /// Evaluate incident velocity (∂p/∂n) at given points with normals
164    ///
165    /// # Arguments
166    /// * `points` - Evaluation points (N × 3 array)
167    /// * `normals` - Unit normal vectors at each point (N × 3 array)
168    /// * `physics` - Physical parameters
169    ///
170    /// # Returns
171    /// Complex normal velocity values (actually ∂p/∂n, not v_n)
172    pub fn evaluate_normal_derivative(
173        &self,
174        points: &Array2<f64>,
175        normals: &Array2<f64>,
176        physics: &PhysicsParams,
177    ) -> Array1<Complex64> {
178        let n = points.nrows();
179        let k = physics.wave_number;
180        let mut dpdn = Array1::zeros(n);
181
182        match self {
183            IncidentField::PlaneWave { direction, amplitude } => {
184                for i in 0..n {
185                    let kdotx = k * (direction[0] * points[[i, 0]]
186                        + direction[1] * points[[i, 1]]
187                        + direction[2] * points[[i, 2]]);
188
189                    // ∂p/∂n = ik (d·n) p
190                    let kdotn = k * (direction[0] * normals[[i, 0]]
191                        + direction[1] * normals[[i, 1]]
192                        + direction[2] * normals[[i, 2]]);
193
194                    let p = *amplitude * Complex64::new(kdotx.cos(), kdotx.sin());
195                    dpdn[i] = Complex64::new(0.0, kdotn) * p;
196                }
197            }
198
199            IncidentField::PointSource { position, strength } => {
200                for i in 0..n {
201                    let dx = points[[i, 0]] - position[0];
202                    let dy = points[[i, 1]] - position[1];
203                    let dz = points[[i, 2]] - position[2];
204                    let r = (dx * dx + dy * dy + dz * dz).sqrt();
205
206                    if r > 1e-10 {
207                        let kr = k * r;
208                        let _r2 = r * r;
209
210                        // ∂G/∂r = (ik - 1/r) * G
211                        let exp_ikr = Complex64::new(kr.cos(), kr.sin());
212                        let g = exp_ikr / (4.0 * PI * r);
213                        let dgdr = (Complex64::new(0.0, k) - Complex64::new(1.0 / r, 0.0)) * g;
214
215                        // ∂r/∂n = (x-x0)·n / r
216                        let drdn = (dx * normals[[i, 0]] + dy * normals[[i, 1]] + dz * normals[[i, 2]]) / r;
217
218                        dpdn[i] = *strength * dgdr * drdn;
219                    }
220                }
221            }
222
223            IncidentField::MultiplePlaneWaves(waves) => {
224                for (direction, amplitude) in waves {
225                    for i in 0..n {
226                        let kdotx = k * (direction[0] * points[[i, 0]]
227                            + direction[1] * points[[i, 1]]
228                            + direction[2] * points[[i, 2]]);
229
230                        let kdotn = k * (direction[0] * normals[[i, 0]]
231                            + direction[1] * normals[[i, 1]]
232                            + direction[2] * normals[[i, 2]]);
233
234                        let p = *amplitude * Complex64::new(kdotx.cos(), kdotx.sin());
235                        dpdn[i] += Complex64::new(0.0, kdotn) * p;
236                    }
237                }
238            }
239
240            IncidentField::MultiplePointSources(sources) => {
241                for (position, strength) in sources {
242                    for i in 0..n {
243                        let dx = points[[i, 0]] - position[0];
244                        let dy = points[[i, 1]] - position[1];
245                        let dz = points[[i, 2]] - position[2];
246                        let r = (dx * dx + dy * dy + dz * dz).sqrt();
247
248                        if r > 1e-10 {
249                            let kr = k * r;
250                            let exp_ikr = Complex64::new(kr.cos(), kr.sin());
251                            let g = exp_ikr / (4.0 * PI * r);
252                            let dgdr = (Complex64::new(0.0, k) - Complex64::new(1.0 / r, 0.0)) * g;
253                            let drdn = (dx * normals[[i, 0]] + dy * normals[[i, 1]] + dz * normals[[i, 2]]) / r;
254
255                            dpdn[i] += *strength * dgdr * drdn;
256                        }
257                    }
258                }
259            }
260        }
261
262        dpdn
263    }
264
265    /// Compute the right-hand side vector for BEM
266    ///
267    /// For exterior Neumann problem (rigid scatterer) using DIRECT formulation
268    /// where the unknown is the total surface pressure p:
269    ///
270    /// The RHS comes from the incident field contribution to the integral equation.
271    /// From NumCalc NC_IncidentWaveRHS:
272    ///   RHS = -(γ*p_inc + β*τ*∂p_inc/∂n)
273    ///
274    /// For exterior problems: τ = +1, γ = 1
275    /// For interior problems: τ = -1
276    pub fn compute_rhs(
277        &self,
278        element_centers: &Array2<f64>,
279        element_normals: &Array2<f64>,
280        physics: &PhysicsParams,
281        use_burton_miller: bool,
282    ) -> Array1<Complex64> {
283        if use_burton_miller {
284            let beta = physics.burton_miller_beta();
285            self.compute_rhs_with_beta(element_centers, element_normals, physics, beta)
286        } else {
287            // CBIE only: RHS = -γ*p_inc
288            let gamma = Complex64::new(physics.gamma(), 0.0);
289            let p_inc = self.evaluate_pressure(element_centers, physics);
290            -gamma * p_inc
291        }
292    }
293
294    /// Compute RHS with custom Burton-Miller coupling parameter
295    ///
296    /// From NumCalc NC_IncidentWaveRHS:
297    ///   RHS = -(γ*p_inc + β*τ*∂p_inc/∂n)
298    ///
299    /// This formula matches the C++ implementation for proper BEM assembly.
300    pub fn compute_rhs_with_beta(
301        &self,
302        element_centers: &Array2<f64>,
303        element_normals: &Array2<f64>,
304        physics: &PhysicsParams,
305        beta: Complex64,
306    ) -> Array1<Complex64> {
307        // Get incident pressure at collocation points
308        let p_inc = self.evaluate_pressure(element_centers, physics);
309
310        // Get incident velocity (normal derivative) at collocation points
311        let dpdn = self.evaluate_normal_derivative(element_centers, element_normals, physics);
312
313        // Parameters from physics
314        let gamma = Complex64::new(physics.gamma(), 0.0);
315        let tau = Complex64::new(physics.tau, 0.0);
316
317        // Burton-Miller RHS: -(γ*p_inc + β*τ*∂p_inc/∂n)
318        // This matches NC_IncidentWaveRHS: zrc -= zg*Gama3 + zv*(zBta3*Tao_)
319        let mut rhs = Array1::zeros(element_centers.nrows());
320        for i in 0..rhs.len() {
321            rhs[i] = -(gamma * p_inc[i] + beta * tau * dpdn[i]);
322        }
323
324        rhs
325    }
326}
327
328#[cfg(test)]
329mod tests {
330    use super::*;
331
332    fn make_physics(k: f64) -> PhysicsParams {
333        let c = 343.0;
334        let freq = k * c / (2.0 * PI);
335        PhysicsParams::new(freq, c, 1.21, false)
336    }
337
338    #[test]
339    fn test_plane_wave_on_axis() {
340        let incident = IncidentField::plane_wave([0.0, 0.0, 1.0], 1.0);
341        let physics = make_physics(1.0);
342
343        // Points on z-axis
344        let points = Array2::from_shape_vec(
345            (3, 3),
346            vec![0.0, 0.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, -1.0],
347        )
348        .unwrap();
349
350        let p = incident.evaluate_pressure(&points, &physics);
351
352        // At z=0: exp(0) = 1
353        assert!((p[0].re - 1.0).abs() < 1e-10);
354        assert!(p[0].im.abs() < 1e-10);
355
356        // At z=1: exp(ik) = exp(i) = cos(1) + i*sin(1)
357        // Since direction is +z, kdotx = k*z = 1 at z=1
358        assert!((p[1].re - (1.0_f64).cos()).abs() < 1e-10);
359        assert!((p[1].im - (1.0_f64).sin()).abs() < 1e-10);
360    }
361
362    #[test]
363    fn test_point_source_decay() {
364        let incident = IncidentField::point_source([0.0, 0.0, 0.0], 1.0);
365        let physics = make_physics(1.0);
366
367        // Points at different distances
368        let points = Array2::from_shape_vec(
369            (3, 3),
370            vec![1.0, 0.0, 0.0, 2.0, 0.0, 0.0, 4.0, 0.0, 0.0],
371        )
372        .unwrap();
373
374        let p = incident.evaluate_pressure(&points, &physics);
375
376        // Pressure should decay as 1/r
377        let ratio_1_2 = p[0].norm() / p[1].norm();
378        let ratio_2_4 = p[1].norm() / p[2].norm();
379
380        assert!((ratio_1_2 - 2.0).abs() < 0.1);
381        assert!((ratio_2_4 - 2.0).abs() < 0.1);
382    }
383
384    #[test]
385    fn test_plane_wave_normal_derivative() {
386        let incident = IncidentField::plane_wave([0.0, 0.0, 1.0], 1.0);
387        let physics = make_physics(1.0);
388
389        // Point with normal pointing in +z direction (along wave direction)
390        let points = Array2::from_shape_vec((1, 3), vec![0.0, 0.0, 0.0]).unwrap();
391        let normals = Array2::from_shape_vec((1, 3), vec![0.0, 0.0, 1.0]).unwrap();
392
393        let dpdn = incident.evaluate_normal_derivative(&points, &normals, &physics);
394
395        // ∂p/∂n = ik (d·n) p = ik * (+1) * 1 = +ik
396        assert!(dpdn[0].re.abs() < 1e-10);
397        assert!((dpdn[0].im - 1.0).abs() < 1e-10);
398    }
399
400    #[test]
401    fn test_rhs_computation() {
402        let incident = IncidentField::plane_wave([0.0, 0.0, 1.0], 1.0);
403        let physics = make_physics(1.0);
404
405        let centers = Array2::from_shape_vec((1, 3), vec![0.0, 0.0, 1.0]).unwrap();
406        let normals = Array2::from_shape_vec((1, 3), vec![0.0, 0.0, 1.0]).unwrap();
407
408        let rhs = incident.compute_rhs(&centers, &normals, &physics, false);
409
410        // RHS = -∂p_inc/∂n, should be non-zero
411        assert!(rhs[0].norm() > 0.0);
412    }
413}