AMDiS  2.10
The Adaptive Multi-Dimensional Simulation Toolbox
VectorBackend.hpp
1 #pragma once
2 
3 #include <boost/numeric/mtl/vector/dense_vector.hpp>
4 
5 #include <dune/common/ftraits.hh>
6 
7 #include <amdis/Output.hpp>
8 #include <amdis/common/FakeContainer.hpp>
9 #include <amdis/linearalgebra/VectorBase.hpp>
10 #include <amdis/typetree/MultiIndex.hpp>
11 
12 namespace AMDiS
13 {
14  class DefaultIndexDistribution;
15 
17  template <class T>
18  class MTLVector
19  : public VectorBase<MTLVector<T>>
20  {
21  public:
23  using BaseVector = mtl::dense_vector<T>;
24 
26  using value_type = typename BaseVector::value_type;
27 
29  using size_type = typename BaseVector::size_type;
30 
31  private:
33  using block_type = value_type;
34 
36  using field_type = typename Dune::FieldTraits<value_type>::field_type;
37 
38  public:
40  explicit MTLVector(DefaultIndexDistribution const&) {}
41 
43  BaseVector const& vector() const
44  {
45  return vector_;
46  }
47 
50  {
51  return vector_;
52  }
53 
55  std::size_t size() const
56  {
57  return mtl::vec::size(vector_);
58  }
59 
60 
62  template <class Basis>
63  void init(Basis const& basis, bool clear)
64  {
65  vector_.change_dim(basis.dimension());
66  if (clear)
67  set_to_zero(vector_);
68  }
69 
71  template <class MultiIndex>
72  value_type& at(MultiIndex const& idx)
73  {
74  const size_type i = flatMultiIndex(idx);
75  test_exit_dbg(i < size(),"Index {} out of range [0,{})", i, size());
76  return vector_[i];
77  }
78 
80  template <class MultiIndex>
81  value_type const& at(MultiIndex const& idx) const
82  {
83  const size_type i = flatMultiIndex(idx);
84  test_exit_dbg(i < size(), "Index {} out of range [0,{})", i, size());
85  return vector_[i];
86  }
87 
88  private:
90  BaseVector vector_;
91  };
92 
93 
94  namespace Recursive
95  {
96  template <class T>
97  struct ForEach<MTLVector<T>> : ForEach<VectorBase<MTLVector<T>>> {};
98 
99  template <class T>
100  struct Transform<MTLVector<T>> : Transform<VectorBase<MTLVector<T>>> {};
101 
102  template <class T>
103  struct InnerProduct<MTLVector<T>> : InnerProduct<VectorBase<MTLVector<T>>> {};
104 
105  } // end namespace Recursive
106 } // end namespace AMDiS
constexpr bool MultiIndex
A multi-index type.
Definition: Concepts.hpp:149
MTLVector(DefaultIndexDistribution const &)
Constructor. Constructs new BaseVector.
Definition: VectorBackend.hpp:40
Definition: ForEach.hpp:17
typename BaseVector::value_type value_type
The type of the elements of the DOFVector.
Definition: VectorBackend.hpp:26
Definition: AdaptBase.hpp:6
Definition: Transform.hpp:15
void init(Basis const &basis, bool clear)
Resize the vector_ to the size s.
Definition: VectorBackend.hpp:63
General implementation of recursive inner-product.
Definition: InnerProduct.hpp:16
value_type & at(MultiIndex const &idx)
Access the entry i of the vector with read-access.
Definition: VectorBackend.hpp:72
mtl::dense_vector< T > BaseVector
The type of the base vector.
Definition: VectorBackend.hpp:23
std::size_t size() const
Return the current size of the vector_.
Definition: VectorBackend.hpp:55
CRTP base class for flat vector backends.
Definition: VectorBase.hpp:18
Dummy implementation for sequential index "distribution".
Definition: IndexDistribution.hpp:6
The basic container that stores a base vector data.
Definition: VectorBackend.hpp:18
BaseVector & vector()
Return the data-vector vector_.
Definition: VectorBackend.hpp:49
value_type const & at(MultiIndex const &idx) const
Access the entry i of the vector with read-access.
Definition: VectorBackend.hpp:81
typename BaseVector::size_type size_type
The index/size - type.
Definition: VectorBackend.hpp:29
BaseVector const & vector() const
Return the data-vector vector_.
Definition: VectorBackend.hpp:43