NiHu  2.0
bbfmm_laplace_3d.cpp
1 #include "library/laplace_3d.hpp"
2 #include "fmm/black_box_fmm.hpp"
3 #include "fmm/cluster_tree.hpp"
6 #include "fmm/fmm_integrated.hpp"
7 #include "fmm/fmm_indexed.hpp"
8 #include "fmm/fmm_precompute.hpp"
9 #include "fmm/fmm_matrix.hpp"
10 #include "fmm/lists.hpp"
12 #include "library/lib_element.hpp"
13 #include "core/mesh.hpp"
14 #include "core/function_space.hpp"
15 #include "core/field.hpp"
16 
25 
28 
29 int main(int argc, char const *argv[])
30 {
31  kernel_t kernel;
32  fmm_t bbfmm(kernel);
33  auto mesh = NiHu::read_off_mesh(argv[1], elem_tag_t());
34  auto const &space = NiHu::constant_view(mesh);
35  cluster_tree_t tree(NiHu::fmm::create_elem_center_iterator(mesh.begin<elem_t>()),
36  NiHu::fmm::create_elem_center_iterator(mesh.end<elem_t>()),
38  NiHu::fmm::interaction_lists lists(tree);
40  field_tag_t(), field_tag_t(), 5, true);
41  auto idx_fctr = NiHu::fmm::create_indexed_functor(
42  space.field_begin<field_t>(),
43  space.field_end<field_t>(),
44  space.field_begin<field_t>(),
45  space.field_end<field_t>(),
46  tree);
47  auto pre_fctr = NiHu::fmm::create_precompute_functor(tree, lists);
48 
49  for (size_t i = 0; i < tree.get_n_clusters(); ++i)
50  tree[i].set_chebyshev_order(5);
51 
53  pre_fctr(idx_fctr(int_fctr(bbfmm.create_p2p()))),
54  pre_fctr(idx_fctr(int_fctr(bbfmm.create_p2m()))),
55  pre_fctr(idx_fctr(int_fctr(bbfmm.create_p2l()))),
56  pre_fctr(idx_fctr(int_fctr(bbfmm.create_m2p()))),
57  pre_fctr(idx_fctr(int_fctr(bbfmm.create_l2p()))),
58  pre_fctr(idx_fctr(bbfmm.create_m2m())),
59  pre_fctr(idx_fctr(bbfmm.create_l2l())),
60  pre_fctr(idx_fctr(bbfmm.create_m2l())),
61  tree,
62  lists);
63 
64  Eigen::Matrix<double, Eigen::Dynamic, 1> xct(mat.cols(), 1);
65  xct.setConstant(1.0);
66  auto res = mat * xct;
67 
68  std::cout << res << std::endl;
69 }
fmm_matrix.hpp
Class NiHu::fmm::fmm_matrix.
NiHu::fmm::black_box_fmm
Black box FMM for a smooth kernel.
Definition: black_box_fmm.hpp:61
NiHu::fmm::chebyshev_cluster
Cluster class of the Black Box FMM.
Definition: chebyshev_cluster.hpp:28
NiHu::fmm::cluster_tree< cluster_t >
read_off_mesh.hpp
export NiHu mesh from OFF format
NiHu::read_off_mesh
mesh< tmp::vector< typename tag2type< Tags >::type... > > read_off_mesh(std::string const &fname, Tags...tags)
Read mesh from OFF format.
Definition: read_off_mesh.hpp:153
NiHu::type2tag
Metafunction assigning a tag to a type.
Definition: type2tag.hpp:17
NiHu::field_view
Field automatically generated from an element using a field generation option.
Definition: field.hpp:222
fmm_operator_collection.hpp
Implementation of NiHu::fmm::fmm_operator_collection.
NiHu::mesh
container class for a mesh
Definition: mesh.hpp:110
NiHu::exponential_covariance_kernel
Definition: covariance_kernel.hpp:42
NiHu::function_space_view
A mesh extended with a Field generating option.
Definition: function_space.hpp:116
field.hpp
implementation of fields and field views
lists.hpp
Definition of class NiHu::fmm::interaction_lists.
NiHu::fmm::divide_num_nodes
Class representing a cluster division based on number of nodes.
Definition: divide.hpp:74
elem_center_iterator.hpp
Classes for iterating through elements and fields.
NiHu::constant_view
const function_space_view< Mesh, field_option::constant, Dimension > & constant_view(mesh_base< Mesh > const &msh, Dimension dim=Dimension())
factory function to transform a mesh into a constant function space
Definition: function_space.hpp:304
NiHu::surface_element
class describing a surface element that provides a normal vector
Definition: element.hpp:451
cluster_tree.hpp
Implementation of class NiHu::fmm::cluster_tree.
fmm_indexed.hpp
Interface for indexing FMM operators.
fmm_precompute.hpp
Operator pre-computation interface.
lib_element.hpp
NiHu::fmm::create_integrated_functor
auto create_integrated_functor(TestTag test_tag, TrialTag trial_tag, size_t quadrature_order, bool sing_check)
factory function to create an integrated functor
Definition: fmm_integrated.hpp:131
black_box_fmm.hpp
Implementation of the Black Box FMM using Chebyshev interpolation.
NiHu::fmm::create_fmm_matrix
auto create_fmm_matrix(fmm_operator_collection< CollOps... > const &collection, cluster_tree< Cluster > const &tree, interaction_lists const &lists)
Factory function to create fmm_matrix from an operator collection.
Definition: fmm_matrix.hpp:660
fmm_integrated.hpp
Integrated operators (P2P, P2X, X2P)
function_space.hpp
declaration of class function_space
NiHu::normal_derivative_kernel
Normal derivative of a distance dependent kernel.
Definition: normal_derivative_kernel.hpp:26
NiHu::fmm::interaction_lists
class storing the different interaction lists of a tree
Definition: lists.hpp:22
mesh.hpp
Declaration of class Mesh.
NiHu::volume_element
class describing a volume element that has no normal vector
Definition: element.hpp:455