NiHu  2.0
helmholtz_2d_wb_x2x_matrix.cpp
Go to the documentation of this file.
1 
8 
9 #include <iostream>
10 
11 namespace NiHu
12 {
13 namespace fmm
14 {
15 
17  : m_level_data_to(nullptr)
18  , m_level_data_from(nullptr)
19 {
20 }
21 
22 
24  helmholtz_2d_wb_level_data const &level_data_to,
25  helmholtz_2d_wb_level_data const &level_data_from,
26  cvector_t diag_coeffs)
27  : m_level_data_to(&level_data_to)
28  , m_level_data_from(&level_data_from)
29  , m_cm(m_level_data_to->get_expansion_length(),
30  m_level_data_from->get_expansion_length(),
31  diag_coeffs)
32 {
33  if (m_level_data_from->is_high_freq() && m_level_data_to->is_high_freq())
34  m_transfer = m_level_data_to->dft(diag_coeffs);
35 }
36 
39 {
40  // low to low
41  if (!m_level_data_from->is_high_freq() && !m_level_data_to->is_high_freq())
42  return m_cm * rhs;
43 
44  // low to high
45  if (!m_level_data_from->is_high_freq() && m_level_data_to->is_high_freq())
46  return m_level_data_to->dft(m_cm * rhs);
47 
48  // high to high
49  return m_level_data_to->interp_up(rhs).array() * m_transfer.array();
50 }
51 
52 
54  : m_level_data_to(nullptr)
55  , m_level_data_from(nullptr)
56 {
57 }
58 
60  helmholtz_2d_wb_level_data const &level_from,
61  cvector_t diag_coeffs)
62  : m_level_data_to(&level_to)
63  , m_level_data_from(&level_from)
64  , m_cm(m_level_data_to->get_expansion_length(),
65  m_level_data_from->get_expansion_length(),
66  diag_coeffs)
67 {
68  if (m_level_data_from->is_high_freq() && m_level_data_to->is_high_freq())
69  m_transfer = m_level_data_from->dft(diag_coeffs);
70 }
71 
74 {
75  if (!m_level_data_from->is_high_freq() && !m_level_data_to->is_high_freq())
76  return m_cm * rhs;
77 
78  // high to low
79  if (m_level_data_from->is_high_freq() && !m_level_data_to->is_high_freq())
80  return m_cm * m_level_data_from->idft(rhs);
81 
82  // L2L high to high
83  return m_level_data_to->interp_dn(m_transfer.array() * rhs.array());
84 }
85 
86 
88  : m_level_data(nullptr)
89 {
90 }
91 
92 
94  cvector_t diag_coeffs)
95  : m_level_data(&level)
96  , m_cm(m_level_data->get_expansion_length(), m_level_data->get_expansion_length(), diag_coeffs)
97 {
98  if (m_level_data->is_high_freq())
99  m_transfer = m_level_data->dft(diag_coeffs);
100 }
101 
104 {
106  if (!m_level_data->is_high_freq())
107  return m_cm * rhs;
108 
110  return rhs.array() * m_transfer.array();
111 }
112 
113 } // end of namespace fmm
114 } // end of namespace NiHu
helmholtz_2d_wb_x2x_matrix.h
Interface of M2M, L2L, and M2L operations for the 2D wide band Helmholtz FMM.
NiHu::fmm::helmholtz_2d_wb_l2l_matrix::helmholtz_2d_wb_l2l_matrix
helmholtz_2d_wb_l2l_matrix()
costructor
Definition: helmholtz_2d_wb_x2x_matrix.cpp:53
NiHu::fmm::helmholtz_2d_wb_level_data::interp_up
const cvector_t & interp_up(cvector_t const &multi) const
interpolate multipole contribution up to this level
Definition: helmholtz_2d_wb_level_data.cpp:166
NiHu::fmm::helmholtz_2d_wb_level_data::dft
const cvector_t & dft(cvector_t const &multi) const
perform dft on this level
Definition: helmholtz_2d_wb_level_data.cpp:120
NiHu::fmm::helmholtz_2d_wb_level_data::is_high_freq
bool is_high_freq() const
return true if the level is in hgh frequency
Definition: helmholtz_2d_wb_level_data.cpp:114
NiHu::fmm::helmholtz_2d_wb_l2l_matrix::cvector_t
Eigen::Matrix< std::complex< double >, Eigen::Dynamic, 1 > cvector_t
complex column vector (the local type)
Definition: helmholtz_2d_wb_x2x_matrix.h:57
NiHu::fmm::helmholtz_2d_wb_level_data
class containing the level data of the helmholtz 2d wide band fmm
Definition: helmholtz_2d_wb_level_data.h:24
NiHu::fmm::helmholtz_2d_wb_m2l_matrix::operator*
cvector_t operator*(cvector_t const &rhs) const
multiply the matrix with a multipole coefficient from the right
Definition: helmholtz_2d_wb_x2x_matrix.cpp:103
NiHu::fmm::helmholtz_2d_wb_m2m_matrix::operator*
cvector_t operator*(cvector_t const &rhs) const
multiply the matrix with a multipole coefficient from the right
Definition: helmholtz_2d_wb_x2x_matrix.cpp:38
NiHu::fmm::helmholtz_2d_wb_level_data::interp_dn
const cvector_t & interp_dn(cvector_t const &local) const
interpolate local data down to this level
Definition: helmholtz_2d_wb_level_data.cpp:180
NiHu::fmm::helmholtz_2d_wb_m2l_matrix::helmholtz_2d_wb_m2l_matrix
helmholtz_2d_wb_m2l_matrix()
costructor
Definition: helmholtz_2d_wb_x2x_matrix.cpp:87
NiHu::fmm::helmholtz_2d_wb_level_data::idft
const cvector_t & idft(cvector_t const &local) const
perform inverse dft on this level
Definition: helmholtz_2d_wb_level_data.cpp:139
NiHu::fmm::helmholtz_2d_wb_m2m_matrix::cvector_t
Eigen::Matrix< std::complex< double >, Eigen::Dynamic, 1 > cvector_t
complex column vector (the multipole type)
Definition: helmholtz_2d_wb_x2x_matrix.h:23
NiHu::fmm::helmholtz_2d_wb_l2l_matrix::operator*
cvector_t operator*(cvector_t const &rhs) const
multiply the matrix with a local coefficient from the right
Definition: helmholtz_2d_wb_x2x_matrix.cpp:73
NiHu::fmm::helmholtz_2d_wb_m2m_matrix::helmholtz_2d_wb_m2m_matrix
helmholtz_2d_wb_m2m_matrix()
costructor
Definition: helmholtz_2d_wb_x2x_matrix.cpp:16
NiHu::fmm::helmholtz_2d_wb_m2l_matrix::cvector_t
Eigen::Matrix< std::complex< double >, Eigen::Dynamic, 1 > cvector_t
complex column vector (the multipole and local type)
Definition: helmholtz_2d_wb_x2x_matrix.h:91