1use ndarray::{Array1, Array2};
7use num_complex::Complex64;
8use std::f64::consts::PI;
9
10use crate::core::integration::gauss::{quad_quadrature, triangle_quadrature};
11use crate::core::integration::singular::generate_subelements;
12use crate::core::mesh::element::{cross_product, normalize};
13use crate::core::types::{ElementType, IntegrationResult, PhysicsParams};
14
15pub fn regular_integration(
34 source_point: &Array1<f64>,
35 source_normal: &Array1<f64>,
36 element_coords: &Array2<f64>,
37 element_type: ElementType,
38 element_area: f64,
39 physics: &PhysicsParams,
40 bc_values: Option<&[Complex64]>,
41 bc_type: i32,
42 compute_rhs: bool,
43) -> IntegrationResult {
44 let wavruim = physics.harmonic_factor * physics.wave_number;
45 let k2 = physics.wave_number * physics.wave_number;
46
47 let mut result = IntegrationResult::default();
48
49 let subelements =
51 generate_subelements(source_point, element_coords, element_type, element_area);
52
53 let (_csi_nodes, _eta_nodes): (Vec<f64>, Vec<f64>) = match element_type {
57 ElementType::Tri3 => (vec![0.0, 1.0, 0.0], vec![0.0, 0.0, 1.0]),
58 ElementType::Quad4 => (vec![1.0, -1.0, -1.0, 1.0], vec![1.0, 1.0, -1.0, -1.0]),
59 };
60
61 for subel in &subelements {
63 let xice = subel.xi_center;
64 let etce = subel.eta_center;
65 let fase = subel.factor;
66 let fase2 = fase * fase;
67 let iforie = (fase.abs() - 1.0).abs() < 1e-10; let quad_points = get_quadrature_points(element_type, subel.gauss_order);
71
72 for (csi_gau, eta_gau, wei_gau) in quad_points {
73 let (xio, eto, weih2) = if iforie {
75 (csi_gau, eta_gau, wei_gau)
76 } else if let Some(tri_verts) = &subel.tri_vertices {
77 let l0 = 1.0 - csi_gau - eta_gau; let xio = tri_verts[0].0 * l0 + tri_verts[1].0 * csi_gau + tri_verts[2].0 * eta_gau;
83 let eto = tri_verts[0].1 * l0 + tri_verts[1].1 * csi_gau + tri_verts[2].1 * eta_gau;
84
85 let dx1 = tri_verts[1].0 - tri_verts[0].0;
88 let dy1 = tri_verts[1].1 - tri_verts[0].1;
89 let dx2 = tri_verts[2].0 - tri_verts[0].0;
90 let dy2 = tri_verts[2].1 - tri_verts[0].1;
91 let det = (dx1 * dy2 - dx2 * dy1).abs();
92 let weih2 = wei_gau * det;
95
96 (xio, eto, weih2)
97 } else {
98 (
100 xice + csi_gau * fase,
101 etce + eta_gau * fase,
102 wei_gau * fase2,
103 )
104 };
105
106 let (shape_fn, jacobian, el_norm, crd_poi) =
108 compute_parameters(element_coords, element_type, xio, eto);
109
110 let wga = weih2 * jacobian;
111
112 let mut diff_fsp = Array1::zeros(3);
114 for i in 0..3 {
115 diff_fsp[i] = crd_poi[i] - source_point[i];
116 }
117
118 let (unit_r, dis_fsp) = normalize(&diff_fsp);
119
120 if dis_fsp < 1e-15 {
121 continue;
122 }
123
124 let re1 = wavruim * dis_fsp;
126 let re2 = wga / (4.0 * PI * dis_fsp);
127 let zg = Complex64::new(re1.cos() * re2, re1.sin() * re2);
128
129 let z1 = Complex64::new(-1.0 / dis_fsp, wavruim);
131 let zhh_base = zg * z1;
132 let re1_h = unit_r.dot(&el_norm); let zhh = zhh_base * re1_h;
134
135 let re2_h = -unit_r.dot(source_normal); let zht = zhh_base * re2_h;
138
139 let rq = re1_h * re2_h;
141 let nx_dot_ny = source_normal.dot(&el_norm);
142 let dq = dis_fsp * dis_fsp;
143
144 let ze_factor = Complex64::new(
145 (3.0 / dq - k2) * rq + nx_dot_ny / dq,
146 -wavruim / dis_fsp * (3.0 * rq + nx_dot_ny),
147 );
148 let ze = zg * ze_factor;
149
150 result.g_integral += zg;
152 result.dg_dn_integral += zhh;
153 result.dg_dnx_integral += zht;
154 result.d2g_dnxdny_integral += ze;
155
156 if compute_rhs && let Some(bc) = bc_values {
158 let mut zbgao = Complex64::new(0.0, 0.0);
160 for (i, &n) in shape_fn.iter().enumerate() {
161 if i < bc.len() {
162 zbgao += bc[i] * n;
163 }
164 }
165
166 let gamma = physics.gamma();
167 let tau = physics.tau;
168 let beta = physics.burton_miller_beta();
169
170 if bc_type == 0 {
171 result.rhs_contribution += (zg * gamma * tau + zht * beta) * zbgao;
173 } else if bc_type == 1 {
174 result.rhs_contribution -= (zhh * gamma * tau + ze * beta) * zbgao;
176 }
177 }
178 }
179 }
180
181 result
182}
183
184fn get_quadrature_points(element_type: ElementType, order: usize) -> Vec<(f64, f64, f64)> {
186 match element_type {
187 ElementType::Tri3 => triangle_quadrature(order),
188 ElementType::Quad4 => quad_quadrature(order),
189 }
190}
191
192fn compute_parameters(
194 element_coords: &Array2<f64>,
195 element_type: ElementType,
196 s: f64,
197 t: f64,
198) -> (Vec<f64>, f64, Array1<f64>, Array1<f64>) {
199 let num_nodes = element_type.num_nodes();
200
201 let (shape_fn, shape_ds, shape_dt) = match element_type {
202 ElementType::Tri3 => {
203 let shape_fn = vec![1.0 - s - t, s, t];
208 let shape_ds = vec![-1.0, 1.0, 0.0];
209 let shape_dt = vec![-1.0, 0.0, 1.0];
210 (shape_fn, shape_ds, shape_dt)
211 }
212 ElementType::Quad4 => {
213 let s1 = 0.25 * (s + 1.0);
215 let s2 = 0.25 * (s - 1.0);
216 let t1 = t + 1.0;
217 let t2 = t - 1.0;
218
219 let shape_fn = vec![s1 * t1, -s2 * t1, s2 * t2, -s1 * t2];
220 let shape_ds = vec![
221 0.25 * (t + 1.0),
222 -0.25 * (t + 1.0),
223 0.25 * (t - 1.0),
224 -0.25 * (t - 1.0),
225 ];
226 let shape_dt = vec![
227 0.25 * (s + 1.0),
228 0.25 * (1.0 - s),
229 0.25 * (s - 1.0),
230 -0.25 * (s + 1.0),
231 ];
232 (shape_fn, shape_ds, shape_dt)
233 }
234 };
235
236 let mut crd_poi = Array1::zeros(3);
238 let mut dx_ds = Array1::zeros(3);
239 let mut dx_dt = Array1::zeros(3);
240
241 for i in 0..num_nodes {
242 for j in 0..3 {
243 crd_poi[j] += shape_fn[i] * element_coords[[i, j]];
244 dx_ds[j] += shape_ds[i] * element_coords[[i, j]];
245 dx_dt[j] += shape_dt[i] * element_coords[[i, j]];
246 }
247 }
248
249 let normal = cross_product(&dx_ds, &dx_dt);
251 let jacobian = normal.dot(&normal).sqrt();
252
253 let el_norm = if jacobian > 1e-15 {
254 normal / jacobian
255 } else {
256 Array1::zeros(3)
257 };
258
259 (shape_fn, jacobian, el_norm, crd_poi)
260}
261
262pub const QUASI_SINGULAR_THRESHOLD: f64 = 3.0;
265
266pub const HIGH_ACCURACY_THRESHOLD: f64 = 2.0;
268
269pub fn optimal_quadrature_order(distance: f64, element_size: f64) -> usize {
276 let ratio = distance / element_size;
277 if ratio < HIGH_ACCURACY_THRESHOLD {
278 4 } else if ratio < QUASI_SINGULAR_THRESHOLD {
280 3 } else {
282 2 }
284}
285
286pub fn quasi_singular_integration(
304 source_point: &Array1<f64>,
305 source_normal: &Array1<f64>,
306 element_coords: &Array2<f64>,
307 element_type: ElementType,
308 element_center: &Array1<f64>,
309 element_area: f64,
310 physics: &PhysicsParams,
311) -> IntegrationResult {
312 let mut dist_sq = 0.0;
314 for i in 0..3 {
315 let d = element_center[i] - source_point[i];
316 dist_sq += d * d;
317 }
318 let distance = dist_sq.sqrt();
319 let element_size = element_area.sqrt();
320 let ratio = distance / element_size;
321
322 if ratio < HIGH_ACCURACY_THRESHOLD {
324 regular_integration(
326 source_point,
327 source_normal,
328 element_coords,
329 element_type,
330 element_area,
331 physics,
332 None,
333 0,
334 false,
335 )
336 } else if ratio < QUASI_SINGULAR_THRESHOLD {
337 let order = optimal_quadrature_order(distance, element_size);
339 regular_integration_fixed_order(
340 source_point,
341 source_normal,
342 element_coords,
343 element_type,
344 physics,
345 order,
346 )
347 } else {
348 regular_integration_fixed_order(
350 source_point,
351 source_normal,
352 element_coords,
353 element_type,
354 physics,
355 2, )
357 }
358}
359
360pub fn regular_integration_fixed_order(
364 source_point: &Array1<f64>,
365 source_normal: &Array1<f64>,
366 element_coords: &Array2<f64>,
367 element_type: ElementType,
368 physics: &PhysicsParams,
369 gauss_order: usize,
370) -> IntegrationResult {
371 let wavruim = physics.harmonic_factor * physics.wave_number;
372 let k2 = physics.wave_number * physics.wave_number;
373
374 let mut result = IntegrationResult::default();
375
376 let quad_points = get_quadrature_points(element_type, gauss_order);
378
379 for (xi, eta, weight) in quad_points {
380 let (_, jacobian, el_norm, crd_poi) =
382 compute_parameters(element_coords, element_type, xi, eta);
383
384 let wga = weight * jacobian;
385
386 let mut diff_fsp = Array1::zeros(3);
388 for i in 0..3 {
389 diff_fsp[i] = crd_poi[i] - source_point[i];
390 }
391
392 let (unit_r, dis_fsp) = normalize(&diff_fsp);
393
394 if dis_fsp < 1e-15 {
395 continue;
396 }
397
398 let kr = wavruim * dis_fsp;
400 let g_scaled = wga / (4.0 * PI * dis_fsp);
401 let zg = Complex64::new(kr.cos() * g_scaled, kr.sin() * g_scaled);
402
403 let z1 = Complex64::new(-1.0 / dis_fsp, wavruim);
405 let zhh_base = zg * z1;
406 let r_dot_ny = unit_r.dot(&el_norm);
407 let r_dot_nx = -unit_r.dot(source_normal);
408 let zhh = zhh_base * r_dot_ny;
409 let zht = zhh_base * r_dot_nx;
410
411 let rq = r_dot_ny * r_dot_nx;
413 let nx_dot_ny = source_normal.dot(&el_norm);
414 let dq = dis_fsp * dis_fsp;
415
416 let ze_factor = Complex64::new(
417 (3.0 / dq - k2) * rq + nx_dot_ny / dq,
418 -wavruim / dis_fsp * (3.0 * rq + nx_dot_ny),
419 );
420 let ze = zg * ze_factor;
421
422 result.g_integral += zg;
424 result.dg_dn_integral += zhh;
425 result.dg_dnx_integral += zht;
426 result.d2g_dnxdny_integral += ze;
427 }
428
429 result
430}
431
432pub fn integrate_g_only(
436 source_point: &Array1<f64>,
437 element_coords: &Array2<f64>,
438 element_type: ElementType,
439 physics: &PhysicsParams,
440 gauss_order: usize,
441) -> Complex64 {
442 let wavruim = physics.harmonic_factor * physics.wave_number;
443 let mut result = Complex64::new(0.0, 0.0);
444
445 let quad_points = get_quadrature_points(element_type, gauss_order);
446
447 for (xi, eta, weight) in quad_points {
448 let (_, jacobian, _, crd_poi) = compute_parameters(element_coords, element_type, xi, eta);
449
450 let wga = weight * jacobian;
451
452 let mut diff = Array1::zeros(3);
453 for i in 0..3 {
454 diff[i] = crd_poi[i] - source_point[i];
455 }
456 let r = diff.dot(&diff).sqrt();
457
458 if r < 1e-15 {
459 continue;
460 }
461
462 let kr = wavruim * r;
463 let g_scaled = wga / (4.0 * PI * r);
464 result += Complex64::new(kr.cos() * g_scaled, kr.sin() * g_scaled);
465 }
466
467 result
468}
469
470pub fn integrate_h_only(
474 source_point: &Array1<f64>,
475 element_coords: &Array2<f64>,
476 element_type: ElementType,
477 physics: &PhysicsParams,
478 gauss_order: usize,
479) -> Complex64 {
480 let wavruim = physics.harmonic_factor * physics.wave_number;
481 let mut result = Complex64::new(0.0, 0.0);
482
483 let quad_points = get_quadrature_points(element_type, gauss_order);
484
485 for (xi, eta, weight) in quad_points {
486 let (_, jacobian, el_norm, crd_poi) =
487 compute_parameters(element_coords, element_type, xi, eta);
488
489 let wga = weight * jacobian;
490
491 let mut diff = Array1::zeros(3);
492 for i in 0..3 {
493 diff[i] = crd_poi[i] - source_point[i];
494 }
495 let r = diff.dot(&diff).sqrt();
496
497 if r < 1e-15 {
498 continue;
499 }
500
501 let kr = wavruim * r;
502 let g_scaled = wga / (4.0 * PI * r);
503 let zg = Complex64::new(kr.cos() * g_scaled, kr.sin() * g_scaled);
504
505 let z1 = Complex64::new(-1.0 / r, wavruim);
506 let r_dot_n: f64 = diff
507 .iter()
508 .zip(el_norm.iter())
509 .map(|(d, n)| d * n)
510 .sum::<f64>()
511 / r;
512
513 result += zg * z1 * r_dot_n;
514 }
515
516 result
517}
518
519#[cfg(test)]
520mod tests {
521 use super::*;
522 use ndarray::{Array2, array};
523
524 fn make_test_triangle() -> Array2<f64> {
525 Array2::from_shape_vec((3, 3), vec![0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 1.0, 0.0]).unwrap()
526 }
527
528 fn make_test_quad() -> Array2<f64> {
529 Array2::from_shape_vec(
530 (4, 3),
531 vec![0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 1.0, 1.0, 0.0, 0.0, 1.0, 0.0],
532 )
533 .unwrap()
534 }
535
536 #[test]
537 fn test_regular_integration_far_field() {
538 let coords = make_test_triangle();
539 let source = array![10.0, 0.0, 0.0]; let normal = array![0.0, 0.0, 1.0];
541
542 let physics = PhysicsParams::new(1000.0, 343.0, 1.21, false);
543
544 let result = regular_integration(
545 &source,
546 &normal,
547 &coords,
548 ElementType::Tri3,
549 0.5,
550 &physics,
551 None,
552 0,
553 false,
554 );
555
556 assert!(result.g_integral.norm().is_finite());
558 assert!(result.g_integral.norm() < 0.1);
559 }
560
561 #[test]
562 fn test_regular_integration_fixed_order() {
563 let coords = make_test_triangle();
564 let source = array![5.0, 5.0, 0.0];
565 let normal = array![0.0, 0.0, 1.0];
566
567 let physics = PhysicsParams::new(1000.0, 343.0, 1.21, false);
568
569 let result = regular_integration_fixed_order(
570 &source,
571 &normal,
572 &coords,
573 ElementType::Tri3,
574 &physics,
575 4,
576 );
577
578 assert!(result.g_integral.norm().is_finite());
579 }
580
581 #[test]
582 fn test_integrate_g_only() {
583 let coords = make_test_triangle();
584 let source = array![2.0, 2.0, 1.0];
585
586 let physics = PhysicsParams::new(1000.0, 343.0, 1.21, false);
587
588 let g = integrate_g_only(&source, &coords, ElementType::Tri3, &physics, 4);
589
590 assert!(g.norm().is_finite());
591 assert!(g.norm() > 0.0);
592 }
593
594 #[test]
595 fn test_quad_element_integration() {
596 let coords = make_test_quad();
597 let source = array![5.0, 0.5, 0.0];
598 let normal = array![0.0, 0.0, 1.0];
599
600 let physics = PhysicsParams::new(1000.0, 343.0, 1.21, false);
601
602 let result = regular_integration_fixed_order(
603 &source,
604 &normal,
605 &coords,
606 ElementType::Quad4,
607 &physics,
608 4,
609 );
610
611 assert!(result.g_integral.norm().is_finite());
612 }
613
614 #[test]
615 fn test_integration_symmetry() {
616 let coords = make_test_triangle();
618 let source1 = array![0.5, 0.5, 1.0];
619 let normal1 = array![0.0, 0.0, 1.0];
620 let source2 = array![0.5, 0.5, -1.0];
621 let normal2 = array![0.0, 0.0, -1.0];
622
623 let physics = PhysicsParams::new(1000.0, 343.0, 1.21, false);
624
625 let result1 = regular_integration_fixed_order(
626 &source1,
627 &normal1,
628 &coords,
629 ElementType::Tri3,
630 &physics,
631 4,
632 );
633
634 let result2 = regular_integration_fixed_order(
635 &source2,
636 &normal2,
637 &coords,
638 ElementType::Tri3,
639 &physics,
640 4,
641 );
642
643 assert!((result1.g_integral.norm() - result2.g_integral.norm()).abs() < 1e-10);
645 }
646
647 #[test]
648 fn test_compute_parameters_triangle() {
649 let coords = make_test_triangle();
650
651 let (shape, jac, normal, pos) = compute_parameters(&coords, ElementType::Tri3, 0.5, 0.25);
652
653 let sum: f64 = shape.iter().sum();
655 assert!((sum - 1.0).abs() < 1e-10);
656
657 assert!((jac - 1.0).abs() < 1e-10);
659
660 assert!(normal[2].abs() > 0.99);
662
663 assert!((pos[0] - 0.5).abs() < 1e-10);
665 assert!((pos[1] - 0.25).abs() < 1e-10);
666 }
667
668 #[test]
669 fn test_compute_parameters_quad() {
670 let coords = make_test_quad();
671
672 let (shape, _jac, normal, _pos) = compute_parameters(&coords, ElementType::Quad4, 0.0, 0.0);
673
674 let sum: f64 = shape.iter().sum();
676 assert!((sum - 1.0).abs() < 1e-10);
677
678 assert!(normal[2].abs() > 0.99);
680 }
681}