7 #ifndef NIHU_HELMHOLTZ_3D_HF_FMM_HPP_INCLUDED
8 #define NIHU_HELMHOLTZ_3D_HF_FMM_HPP_INCLUDED
15 #include "helmholtz_3d_hf_shift.h"
21 #include "library/helmholtz_3d.hpp"
23 #include <boost/math/constants/constants.hpp>
24 #include <boost/math/special_functions/hankel.hpp>
25 #include <boost/math/special_functions/legendre.hpp>
27 #include <Eigen/Dense>
37 template <
class WaveNumber>
44 using cvector_t = Eigen::Matrix<std::complex<double>, Eigen::Dynamic, 1>;
69 using namespace boost::math::double_constants;
70 auto kd = two_pi * drel;
71 return size_t(std::ceil(kd +
C * std::log(kd + pi)));
79 m_accuracy = accuracy;
86 using namespace boost::math::double_constants;
87 double lambda = two_pi / std::real(m_wave_number);
88 m_level_data_vector.clear();
94 auto &ld = m_level_data_vector[i];
97 auto d = tree[idx].get_bounding_box().get_diameter();
99 ld.set_expansion_length(L);
105 auto &ld = m_level_data_vector[i];
106 auto const &Sto = ld.get_unit_sphere();
108 ld.set_interp_dn(
interpolator(m_level_data_vector[i - 1].get_unit_sphere(), Sto));
110 ld.set_interp_up(
interpolator(m_level_data_vector[i + 1].get_unit_sphere(), Sto));
119 return m_level_data_vector[idx];
127 return m_level_data_vector[idx];
140 m2m(wave_number_t
const &wave_number)
147 return bounding_box_t::dist2idx(
154 using namespace std::complex_literals;
158 auto const &k = this->get_wave_number();
159 cvector_t shift = exp((-1.i * k * Sto.get_s().transpose() * D).array());
175 l2l(wave_number_t
const &wave_number)
182 return bounding_box_t::dist2idx(
189 using namespace std::complex_literals;
197 auto const &k = this->get_wave_number();
199 cvector_t shift = exp((-1.i * k * Sfrom.get_s().transpose() * D).array());
201 cvector_t shift = exp((-1.i * k * Sto.get_s().transpose() * D).array());
224 p2m(wave_number_t
const &wave_number)
234 result_t operator()(
test_input_t const &to, trial_input_t
const &tri)
const
236 return eval(to, tri, std::integral_constant<int, Ny>());
241 result_t eval(
test_input_t const &to, trial_input_t
const &tri,
242 std::integral_constant<int, 0>)
const
244 using namespace std::complex_literals;
246 auto const &y = tri.get_x();
248 auto const &k = this->get_wave_number();
250 return Eigen::exp(-1.i * k * (s.transpose() * d).array());
254 result_t eval(
test_input_t const &to, trial_input_t
const &tri,
255 std::integral_constant<int, 1>)
const
257 using namespace std::complex_literals;
259 auto const &y = tri.get_x();
261 auto const &k = this->get_wave_number();
263 return Eigen::exp(-1.i * k * (s.transpose() * d).array())
264 * ((1.i * k) * (s.transpose() * tri.get_unit_normal()).array());
284 p2l(wave_number_t
const &wave_number)
294 result_t operator()(
test_input_t const &to, trial_input_t
const &tri)
const
296 return eval(to, tri, std::integral_constant<int, Ny>());
300 result_t eval(
test_input_t const &to, trial_input_t
const &tri,
301 std::integral_constant<int, 0>)
const
303 throw std::logic_error(
"Unimplemented p2l operator");
306 result_t eval(
test_input_t const &to, trial_input_t
const &tri,
307 std::integral_constant<int, 1>)
const
309 throw std::logic_error(
"Unimplemented p2l operator");
328 using result_t = Eigen::Matrix<std::complex<double>, 1, Eigen::Dynamic>;
330 l2p(wave_number_t
const &wave_number)
340 result_t operator()(test_input_t
const &tsi,
trial_input_t const &from)
const
342 return eval(tsi, from, std::integral_constant<int, Nx>());
346 result_t eval(test_input_t
const &tsi,
trial_input_t const &from,
347 std::integral_constant<int, 0>)
const
349 using namespace std::complex_literals;
351 auto const &x = tsi.get_x();
353 auto const &s = S.
get_s();
354 auto const &w = S.get_w();
355 auto const &k = this->get_wave_number();
357 return Eigen::exp(-1.i * k * (s.transpose() * d).array()) * w.array();
360 result_t eval(test_input_t
const &tsi,
trial_input_t const &from,
361 std::integral_constant<int, 1>)
const
363 using namespace std::complex_literals;
365 auto const &x = tsi.get_x();
367 auto const &s = S.
get_s();
368 auto const &w = S.get_w();
369 auto const &k = this->get_wave_number();
371 return Eigen::exp(-1.i * k * (s.transpose() * d).array())
372 * (-1.i * k) * (s.transpose() * tsi.get_unit_normal()).array()
392 using result_t = Eigen::Matrix<std::complex<double>, 1, Eigen::Dynamic>;
394 m2p(wave_number_t
const &wave_number)
404 result_t operator()(test_input_t
const &tsi,
trial_input_t const &from)
const
406 return eval(tsi, from, std::integral_constant<int, Nx>());
410 result_t eval(test_input_t
const &tsi,
trial_input_t const &from,
411 std::integral_constant<int, 0>)
const
413 throw std::logic_error(
"Unimplemented m2p operator");
416 result_t eval(test_input_t
const &tsi,
trial_input_t const &from,
417 std::integral_constant<int, 1>)
const
419 throw std::logic_error(
"Unimplemented m2p operator");
430 using result_t = Eigen::DiagonalMatrix<std::complex<double>, Eigen::Dynamic>;
445 auto const &k = this->get_wave_number();
454 return m2l_matrix_impl(Dvec, s, k, L).asDiagonal();
459 Eigen::Matrix<double, 3, Eigen::Dynamic>
const &s,
463 using namespace boost::math::double_constants;
464 using namespace std::complex_literals;
465 using boost::math::legendre_p;
466 using boost::math::sph_hankel_1;
468 double D = Dvec.norm();
472 Eigen::Matrix<double, 1, Eigen::Dynamic> x = Uvec.transpose() * s;
477 for (
size_t l = 0; l <= L; ++l)
479 auto h = sph_hankel_1(l, z);
480 auto c = double(2 * l + 1) * std::pow(1.i, l) * h;
481 for (
int i = 0; i < s.cols(); ++i)
482 M2L(i) += c * legendre_p(
int(l), x(0, i));
484 M2L *= -1.i * k / (4.0 * pi) / (4.0 * pi);
513 template <
int Nx,
int Ny>
560 template <
int Nx,
int Ny>
573 return m2m(m_wave_number);
580 return l2l(m_wave_number);
587 return m2l(m_wave_number);
591 wave_number_t m_wave_number;
592 std::vector<helmholtz_3d_hf_level_data> m_level_data_vector;