dune-functions 2.9.0
flatvectorview.hh
Go to the documentation of this file.
1// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=4 sw=2 sts=2:
3#ifndef DUNE_FUNCTIONS_FUNCTIONSPACEBASES_FLATVECTORVIEW_HH
4#define DUNE_FUNCTIONS_FUNCTIONSPACEBASES_FLATVECTORVIEW_HH
5
6
7#include <array>
8
9#include <dune/common/concept.hh>
10#include <dune/common/hybridutilities.hh>
11#include <dune/common/indices.hh>
12
14
15
16
17
18namespace Dune {
19namespace Functions {
20namespace Impl {
21
22
23template<class V>
24struct FlatVectorBackend
25{
26
27 template<class VV, class Index,
28 typename std::enable_if< models<Concept::HasIndexAccess, VV, Index>(), int>::type = 0>
29 static decltype(auto) getEntry(VV&& v, const Index& i)
30 {
31 return v[i];
32 }
33
34 template<class VV, class Index,
35 typename std::enable_if< not models<Concept::HasIndexAccess, VV, Index>(), int>::type = 0>
36 static decltype(auto) getEntry(VV&& v, const Index&)
37 {
38 return std::forward<VV>(v);
39 }
40
41 template<class VV,
42 typename std::enable_if< models<Concept::HasSizeMethod, VV>(), int>::type = 0>
43 static auto size(VV&& v)
44 {
45 return Dune::Hybrid::size(v);
46 }
47
48 template<class VV,
49 typename std::enable_if< not models<Concept::HasSizeMethod, VV>(), int>::type = 0>
50 static auto size(VV&&)
51 {
52 return Dune::index_constant<1>{};
53 }
54};
55
56
57
58
59template<class K, int n, int m>
60struct FlatVectorBackend<typename Dune::FieldMatrix<K, n, m> >
61{
62
63 template<class VV, class Index>
64 static decltype(auto) getEntry(VV&& v, const Index& i)
65 {
66 return v[i/m][i%m];
67 }
68
69 template<class VV>
70 static auto size(VV&& v)
71 {
72 return Dune::index_constant<n*m>{};
73 }
74};
75
76
77
78template<class K, std::size_t n>
79struct FlatVectorBackend< std::array<K, n> >
80{
81
82 template<class VV, class Index>
83 static decltype(auto) getEntry(VV&& v, const Index& i)
84 {
85 const auto innerSize = decltype(FlatVectorBackend<K>::size(v[0]))::value;
86 return FlatVectorBackend<K>::getEntry(v[i/innerSize], i%innerSize);
87 }
88
89 template<class VV>
90 static auto size(VV&& v)
91 {
92 const auto innerSize = decltype(FlatVectorBackend<K>::size(v[0]))::value;
93 return Dune::index_constant<n*innerSize>{};
94 }
95
96};
97
98
99
100
101template<class T>
102class FlatVectorView
103{
104 using Backend = FlatVectorBackend<std::decay_t<T>>;
105public:
106 FlatVectorView(T& t) :
107 t_(&t)
108 {}
109
110 auto size() const
111 {
112 return Backend::size(*t_);
113 }
114
115 template<class Index>
116 decltype(auto) operator[](const Index& i) const
117 {
118 return Backend::getEntry(*t_, i);
119 }
120
121 template<class Index>
122 decltype(auto) operator[](const Index& i)
123 {
124 return Backend::getEntry(*t_, i);
125 }
126
127private:
128 T* t_;
129};
130
131
132template<class T>
133class FlatVectorView<T&&>
134{
135 using Backend = FlatVectorBackend<std::decay_t<T>>;
136public:
137 FlatVectorView(T&& t) :
138 t_(std::move(t))
139 {}
140
141 auto size() const
142 {
143 return Backend::size(t_);
144 }
145
146 template<class Index>
147 decltype(auto) operator[](const Index& i) const
148 {
149 return Backend::getEntry(t_, i);
150 }
151
152 template<class Index>
153 decltype(auto) operator[](const Index& i)
154 {
155 return Backend::getEntry(t_, i);
156 }
157
158private:
159 T t_;
160};
161
162} // namespace Impl
163
164
165
178template<class T>
180{
181 return Impl::FlatVectorView<T>(t);
182}
183
196template<class T>
197auto flatVectorView(const T& t)
198{
199 return Impl::FlatVectorView<const T>(t);
200}
201
214template<class T>
215auto flatVectorView(T&& t)
216{
217 return Impl::FlatVectorView<T&&>(std::move(t));
218}
219
220
221} // namespace Dune::Functions
222} // namespace Dune
223
224
225#endif // DUNE_FUNCTIONS_FUNCTIONSPACEBASES_FLATVECTORVIEW_HH
Definition: polynomial.hh:10
auto flatVectorView(T &t)
Create flat vector view of passed mutable container.
Definition: flatvectorview.hh:179