NiHu  2.0
helmholtz_2d_wb_level_data.h
Go to the documentation of this file.
1 
7 #ifndef FMM_HELMHOLTZ_2D_WB_LEVEL_DATA_HPP_INCLUDED
8 #define FMM_HELMHOLTZ_2D_WB_LEVEL_DATA_HPP_INCLUDED
9 
10 #include "spectral_interpolate.h"
11 
12 #include <fftw3.h>
13 #include <Eigen/Dense>
14 
15 #include <complex>
16 #include <vector>
17 
18 namespace NiHu
19 {
20 namespace fmm
21 {
22 
25 {
26 public:
28  typedef Eigen::Matrix<std::complex<double>, Eigen::Dynamic, 1> cvector_t;
29 
38 
42  void init(double drel);
43 
47  void set_expansion_length(size_t expansion_length);
48 
50  size_t get_expansion_length() const;
51 
53  size_t get_data_size() const;
54 
56  void set_high_freq(bool hf);
57 
59  bool is_high_freq() const;
60 
64  void set_interp_up(size_t Nfrom);
65 
69  void set_interp_dn(size_t Nfrom);
70 
75  cvector_t const &interp_up(cvector_t const &multi) const;
76 
81  cvector_t const &interp_dn(cvector_t const &local) const;
82 
87  cvector_t const &dft(cvector_t const &multi) const;
88 
93  cvector_t const &idft(cvector_t const &local) const;
94 
95  void set_num_threads(size_t n);
96 
97 private:
98  size_t m_expansion_length;
99  bool m_high_freq;
100 
101  std::vector<spectral_interpolate> m_interp_ups; // should be thread private
102  std::vector<spectral_interpolate> m_interp_dns; // should be thread private
103 
104  fftw_plan m_dft_plan;
105  mutable std::vector<cvector_t> m_multi_paddeds; // should be thread private
106  mutable std::vector<cvector_t> m_multi_spectrums; // should be thread private
107 
108  fftw_plan m_idft_plan;
109  mutable std::vector <cvector_t> m_local_paddeds; // should be thread private
110  mutable std::vector <cvector_t> m_locals; // should be thread private
111 };
112 
113 } // end of namespace fmm
114 } // end of namespace NiHu
115 
116 #endif /* FMM_HELMHOLTZ_2D_WB_LEVEL_DATA_HPP_INCLUDED */
NiHu::fmm::helmholtz_2d_wb_level_data::init
void init(double drel)
initialize the data for a given frequency
Definition: helmholtz_2d_wb_level_data.cpp:52
NiHu::fmm::helmholtz_2d_wb_level_data::get_data_size
size_t get_data_size() const
return the size of the multipole / local data
Definition: helmholtz_2d_wb_level_data.cpp:100
NiHu::fmm::helmholtz_2d_wb_level_data::cvector_t
Eigen::Matrix< std::complex< double >, Eigen::Dynamic, 1 > cvector_t
complex column vector type
Definition: helmholtz_2d_wb_level_data.h:28
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
spectral_interpolate.h
Interface of class NiHu::fmm::spectral_interpolate.
NiHu::fmm::helmholtz_2d_wb_level_data::set_interp_dn
void set_interp_dn(size_t Nfrom)
compute the downward interpolating matrix
Definition: helmholtz_2d_wb_level_data.cpp:172
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_level_data::set_expansion_length
void set_expansion_length(size_t expansion_length)
set the expansion length of the level
Definition: helmholtz_2d_wb_level_data.cpp:60
NiHu::fmm::helmholtz_2d_wb_level_data::set_high_freq
void set_high_freq(bool hf)
set if the level is in high frequency
Definition: helmholtz_2d_wb_level_data.cpp:109
NiHu::fmm::helmholtz_2d_wb_level_data::set_interp_up
void set_interp_up(size_t Nfrom)
compute upward interpolating matrix
Definition: helmholtz_2d_wb_level_data.cpp:158
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_level_data::helmholtz_2d_wb_level_data
helmholtz_2d_wb_level_data()
constructor
Definition: helmholtz_2d_wb_level_data.cpp:16
NiHu::fmm::helmholtz_2d_wb_level_data::~helmholtz_2d_wb_level_data
~helmholtz_2d_wb_level_data()
destructor
Definition: helmholtz_2d_wb_level_data.cpp:29
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_level_data::get_expansion_length
size_t get_expansion_length() const
return the expansion length
Definition: helmholtz_2d_wb_level_data.cpp:95
NiHu::fmm::helmholtz_2d_wb_level_data::operator=
helmholtz_2d_wb_level_data & operator=(helmholtz_2d_wb_level_data const &)=delete
assignment operator