NiHu
2.0
|
the fmm for the 3D Helmholtz equation More...
#include <helmholtz_3d_hf_fmm.hpp>
Classes | |
class | l2l |
L2L operator of the FMM for the Helmholtz equation in 3D. More... | |
class | l2p |
L2P operator of the FMM for the Helmholtz equation in 3D. More... | |
class | m2l |
M2L operator of the FMM for the Helmholtz equation in 3D. More... | |
class | m2m |
M2M operator of the FMM for the Helmholtz equation in 3D. More... | |
class | m2p |
M2P operator of the FMM for the Helmholtz equation in 3D. More... | |
class | p2l |
P2L operator of the FMM for the Helmholtz equation in 3D. More... | |
class | p2m |
P2M operator of the FMM for the Helmholtz equation in 3D. More... | |
Public Types | |
using | wave_number_t = WaveNumber |
template argument as nested type | |
using | cvector_t = Eigen::Matrix< std::complex< double >, Eigen::Dynamic, 1 > |
complex dynamic vector | |
using | cluster_t = helmholtz_3d_hf_cluster |
the fmm's cluster type | |
using | bounding_box_t = typename cluster_t::bounding_box_t |
the bounding box type | |
using | location_t = typename bounding_box_t::location_t |
the physical location type | |
using | cluster_tree_t = cluster_tree< cluster_t > |
the cluster tree type | |
using | distance_dependent_kernel_t = NiHu::helmholtz_kernel< NiHu::space_3d<>, wave_number_t > |
template<int Ny> | |
using | p2m_type = p2m< Ny > |
metafunction returning the P2M operator's type for a given order More... | |
template<int Ny> | |
using | p2l_type = p2l< Ny > |
metafunction returning the P2L operator's type for a given order More... | |
template<int Nx> | |
using | m2p_type = m2p< Nx > |
metafunction returning the M2P operator's type for a given order More... | |
template<int Nx> | |
using | l2p_type = l2p< Nx > |
metafunction returning the L2P operator's type for a given order More... | |
template<int Nx, int Ny> | |
using | p2p_type = fmm::p2p< NiHu::normal_derivative_kernel< distance_dependent_kernel_t, Nx, Ny > > |
metafunction returning the P2P operator's type More... | |
Public Member Functions | |
helmholtz_3d_hf_fmm (wave_number_t const &k) | |
constructor More... | |
void | set_accuracy (double accuracy) |
set the method's accuracy parameter More... | |
void | init_level_data (cluster_tree_t const &tree) |
initialize the level data of the fmm method More... | |
const helmholtz_3d_hf_level_data & | get_level_data (size_t idx) const |
return level data for a specific level More... | |
helmholtz_3d_hf_level_data & | get_level_data (size_t idx) |
return level data reference for a specific level More... | |
template<int Ny> | |
p2m_type< Ny > | create_p2m () const |
factory function for the P2M operator More... | |
template<int Ny> | |
p2l_type< Ny > | create_p2l () const |
factory function for the P2L operator More... | |
template<int Nx> | |
l2p_type< Nx > | create_l2p () const |
factory function for the L2P operator More... | |
template<int Nx> | |
m2p_type< Nx > | create_m2p () const |
factory function for the M2P operator More... | |
template<int Nx, int Ny> | |
p2p_type< Nx, Ny > | create_p2p () const |
factory function for the P2P operator More... | |
m2m | create_m2m () const |
factory function for the M2M operator More... | |
l2l | create_l2l () const |
factory function for the L2L operator More... | |
m2l | create_m2l () const |
factory function for the M2L operator More... | |
Static Public Member Functions | |
static size_t | compute_expansion_length (double drel, double C) |
compute expansion length based on relative cluster size More... | |
the fmm for the 3D Helmholtz equation
WaveNumber | the wave number type (double or complex) |
Definition at line 101 of file helmholtz_3d_hf_fmm.hpp.
using NiHu::fmm::helmholtz_3d_hf_fmm< WaveNumber >::l2p_type = l2p<Nx> |
metafunction returning the L2P operator's type for a given order
Nx | the order of receiver differentiation |
Definition at line 569 of file helmholtz_3d_hf_fmm.hpp.
using NiHu::fmm::helmholtz_3d_hf_fmm< WaveNumber >::m2p_type = m2p<Nx> |
metafunction returning the M2P operator's type for a given order
Nx | the order of receiver differentiation |
Definition at line 564 of file helmholtz_3d_hf_fmm.hpp.
using NiHu::fmm::helmholtz_3d_hf_fmm< WaveNumber >::p2l_type = p2l<Ny> |
metafunction returning the P2L operator's type for a given order
Ny | the order of source differentiation |
Definition at line 559 of file helmholtz_3d_hf_fmm.hpp.
using NiHu::fmm::helmholtz_3d_hf_fmm< WaveNumber >::p2m_type = p2m<Ny> |
metafunction returning the P2M operator's type for a given order
Ny | the order of source differentiation |
Definition at line 554 of file helmholtz_3d_hf_fmm.hpp.
using NiHu::fmm::helmholtz_3d_hf_fmm< WaveNumber >::p2p_type = fmm::p2p< NiHu::normal_derivative_kernel< distance_dependent_kernel_t, Nx, Ny > > |
metafunction returning the P2P operator's type
Nx | the order of receiver differentiation |
Ny | the order of source differentiation |
Definition at line 579 of file helmholtz_3d_hf_fmm.hpp.
|
inline |
|
inlinestatic |
compute expansion length based on relative cluster size
[in] | drel | relative cluster diameter |
[in] | C | accuracy parameter |
Definition at line 130 of file helmholtz_3d_hf_fmm.hpp.
|
inline |
factory function for the L2L operator
Definition at line 639 of file helmholtz_3d_hf_fmm.hpp.
|
inline |
factory function for the L2P operator
Nx | the order of receiver differentiation |
Definition at line 603 of file helmholtz_3d_hf_fmm.hpp.
|
inline |
factory function for the M2L operator
Definition at line 646 of file helmholtz_3d_hf_fmm.hpp.
|
inline |
factory function for the M2M operator
Definition at line 632 of file helmholtz_3d_hf_fmm.hpp.
|
inline |
factory function for the M2P operator
Nx | the order of receiver differentiation |
Definition at line 612 of file helmholtz_3d_hf_fmm.hpp.
|
inline |
factory function for the P2L operator
Ny | the order of source differentiation |
Definition at line 594 of file helmholtz_3d_hf_fmm.hpp.
|
inline |
factory function for the P2M operator
Ny | the order of source differentiation |
Definition at line 585 of file helmholtz_3d_hf_fmm.hpp.
|
inline |
factory function for the P2P operator
Ny | the order of source differentiation |
Nx | the order of receiver differentiation |
Definition at line 622 of file helmholtz_3d_hf_fmm.hpp.
|
inline |
return level data reference for a specific level
[in] | idx | the level index |
Definition at line 188 of file helmholtz_3d_hf_fmm.hpp.
|
inline |
return level data for a specific level
[in] | idx | the level index |
Definition at line 180 of file helmholtz_3d_hf_fmm.hpp.
|
inline |
initialize the level data of the fmm method
[in] | tree | the cluster tree |
Definition at line 147 of file helmholtz_3d_hf_fmm.hpp.
|
inline |
set the method's accuracy parameter
[in] | C | the accuracy parameter |
the accuracy parameter is usually set to 3.0
Definition at line 140 of file helmholtz_3d_hf_fmm.hpp.