NiHu  2.0
p2x_cluster_indexed.hpp
Go to the documentation of this file.
1 
7 #ifndef NIHU_P2X_CLUSTER_INDEXED_HPP_INCLUDED
8 #define NIHU_P2X_CLUSTER_INDEXED_HPP_INCLUDED
9 
10 #include "cluster_tree.hpp"
11 #include "fmm_operator.hpp"
12 #include "util/matrix_traits.hpp"
13 
14 #include <type_traits>
15 
16 namespace NiHu
17 {
18 namespace fmm
19 {
20 
21 template <class Operator>
23  : public fmm_operator<typename std::decay<Operator>::type::fmm_tag>
24 {
25 public:
26  typedef typename std::decay<Operator>::type operator_t;
27  typedef typename operator_t::test_input_t cluster_t;
29  typedef typename operator_t::result_t op_result_t;
30  typedef typename scalar<op_result_t>::type scalar_t;
31  typedef Eigen::Matrix<scalar_t, Eigen::Dynamic, Eigen::Dynamic> result_t;
32  static size_t const cols = num_cols<op_result_t>::value;
33  static size_t const num_dof_per_src = cols;
34 
35  p2x_cluster_indexed(Operator &&op, tree_t const &tree)
36  : m_op(std::forward<Operator>(op))
37  , m_tree(tree)
38  {
39  }
40 
41  result_t operator()(size_t to, size_t from) const
42  {
44  bool first = true;
45  size_t rows = 0;
46  cluster_t const &clus_to = m_tree[to];
47  cluster_t const &clus_from = m_tree[from];
48 
49  result_t mat;
50  for (size_t jj = 0; jj < clus_from.get_n_src_nodes(); ++jj)
51  {
52  size_t j = clus_from.get_src_node_idx()[jj];
53  typename operator_t::result_t res = m_op(clus_to, j);
54  if (first)
55  {
56  rows = res.rows();
57  mat.resize(rows, clus_from.get_n_src_nodes() * cols);
58  first = false;
59  }
60  mat.block(0, cols * jj, rows, cols) = res;
61  }
62 
63  return mat;
64  }
65 
66  result_t operator()(size_t idx) const
67  {
68  return (*this)(idx, idx);
69  }
70 
71  tree_t const &get_tree() const
72  {
73  return m_tree;
74  }
75 
76 private:
77  Operator m_op;
78  tree_t const &m_tree;
79 };
80 
81 
82 template <class Operator>
83 p2x_cluster_indexed<Operator>
84 create_p2x_cluster_indexed(Operator &&op,
85  cluster_tree<typename std::decay<Operator>::type::test_input_t> const &tree)
86 {
87  return p2x_cluster_indexed<Operator>(std::forward<Operator>(op), tree);
88 }
89 
90 } // end of namespace fmm
91 } // namespace NiHu
92 
93 #endif /* NIHU_P2X_CLUSTER_INDEXED_HPP_INCLUDED */
NiHu::fmm::cluster_tree< cluster_t >
NiHu::num_cols
metafunction returning the number of compile time columns
Definition: matrix_traits.hpp:41
matrix_traits.hpp
compile time properties of matrices
cluster_tree.hpp
Implementation of class NiHu::fmm::cluster_tree.
fmm_operator.hpp
FMM operator types and tags.
NiHu::fmm::p2x_cluster_indexed::operator()
result_t operator()(size_t to, size_t from) const
Definition: p2x_cluster_indexed.hpp:41
NiHu::fmm::fmm_operator
Operator defining its tag type.
Definition: fmm_operator.hpp:85
NiHu::fmm::p2x_cluster_indexed
Definition: p2x_cluster_indexed.hpp:22