DUNE PDELab (git)

uncachedmatrixview.hh
1// -*- tab-width: 4; indent-tabs-mode: nil -*-
2#ifndef DUNE_PDELAB_BACKEND_COMMON_UNCACHEDMATRIXVIEW_HH
3#define DUNE_PDELAB_BACKEND_COMMON_UNCACHEDMATRIXVIEW_HH
4
5#include <type_traits>
6
7namespace Dune {
8 namespace PDELab {
9
10
11 template<typename M_, typename RowCache, typename ColCache>
12 class ConstUncachedMatrixView
13 {
14
15 public:
16
17 typedef typename std::remove_const<M_>::type Container;
18
19 static_assert(
20 (std::is_same<
21 typename RowCache::LocalFunctionSpace::Traits::GridFunctionSpace,
22 typename Container::TestGridFunctionSpace
23 >::value),
24 "The RowCache passed to LocalView must belong to the underlying GFSV"
25 );
26
27 static_assert(
28 (std::is_same<
29 typename ColCache::LocalFunctionSpace::Traits::GridFunctionSpace,
30 typename Container::TrialGridFunctionSpace
31 >::value),
32 "The ColCache passed to LocalView must belong to the underlying GFSU"
33 );
34
35 public:
36
37 typedef typename Container::field_type E;
38 typedef typename Container::size_type size_type;
39
40 typedef E ElementType;
41
42 typedef RowCache RowIndexCache;
43 typedef ColCache ColIndexCache;
44
45 typedef typename RowCache::LocalFunctionSpace LFSV;
46 typedef typename ColCache::LocalFunctionSpace LFSU;
47
48 typedef typename LFSV::Traits::DOFIndex RowDOFIndex;
49 typedef typename LFSV::Traits::GridFunctionSpace::Ordering::Traits::ContainerIndex RowContainerIndex;
50
51 typedef typename LFSU::Traits::DOFIndex ColDOFIndex;
52 typedef typename LFSU::Traits::GridFunctionSpace::Ordering::Traits::ContainerIndex ColContainerIndex;
53
54 ConstUncachedMatrixView()
55 : _container(nullptr)
56 , _row_cache(nullptr)
57 , _col_cache(nullptr)
58 {}
59
60 ConstUncachedMatrixView(M_& container)
61 : _container(&container)
62 , _row_cache(nullptr)
63 , _col_cache(nullptr)
64 {}
65
66 const RowIndexCache& rowIndexCache() const
67 {
68 assert(_row_cache);
69 return *_row_cache;
70 }
71
72 const ColIndexCache& colIndexCache() const
73 {
74 assert(_col_cache);
75 return *_col_cache;
76 }
77
78 void attach(M_& container)
79 {
80 _container = &container;
81 }
82
83 void detach()
84 {
85 _container = nullptr;
86 }
87
88 void bind(const RowCache& row_cache, const ColCache& col_cache)
89 {
90 _row_cache = &row_cache;
91 _col_cache = &col_cache;
92 }
93
94 void unbind()
95 {}
96
97 size_type N() const
98 {
99 return rowIndexCache().size();
100 }
101
102 size_type M() const
103 {
104 return colIndexCache().size();
105 }
106
107 template<typename LC>
108 void read(LC& local_container) const
109 {
110 for (size_type i = 0; i < N(); ++i)
111 for (size_type j = 0; j < M(); ++j)
112 local_container.getEntry(i,j) = container()(rowIndexCache().containerIndex(i),colIndexCache().containerIndex(j));
113 }
114
115
116
117 const ElementType& operator()(size_type i, size_type j) const
118 {
119 return container()(rowIndexCache().containerIndex(i),colIndexCache().containerIndex(j));
120 }
121
122 // disable this function if DOFIndex and ContainerIndex have the same type - required for interoperability
123 // with function spaces based on dune-functions bases
124 template<typename RDI, typename CDI>
125 std::enable_if_t<
126 (std::is_same<RDI,RowDOFIndex>{} and std::is_same<CDI,ColDOFIndex>{} and not
127 (std::is_same<RDI,RowContainerIndex>{} and std::is_same<CDI,ColContainerIndex>{})),
128 const ElementType&
129 >
130 operator()(const RDI& i, const CDI& j) const
131 {
132 return container()(rowIndexCache().containerIndex(i),colIndexCache().containerIndex(j));
133 }
134
135 const ElementType& operator()(const RowContainerIndex& i, const ColContainerIndex& j) const
136 {
137 return container()(i,j);
138 }
139
140 const ElementType& operator()(const RowContainerIndex& i, size_type j) const
141 {
142 return container()(i,colIndexCache().containerIndex(j));
143 }
144
145 const ElementType& operator()(size_type i, const ColContainerIndex& j) const
146 {
147 return container()(rowIndexCache().containerIndex(i),j);
148 }
149
150 const Container& container() const
151 {
152 return *_container;
153 }
154
155 protected:
156
157 M_* _container;
158 const RowCache* _row_cache;
159 const ColCache* _col_cache;
160
161 };
162
163
164 template<typename M_, typename RowCache, typename ColCache>
165 class UncachedMatrixView
166 : public ConstUncachedMatrixView<M_,RowCache,ColCache>
167 {
168
169 typedef ConstUncachedMatrixView<M_,RowCache,ColCache> BaseT;
170
171 public:
172
173 typedef M_ Container;
174 typedef typename Container::ElementType ElementType;
175 typedef typename Container::size_type size_type;
176
177 typedef RowCache RowIndexCache;
178 typedef ColCache ColIndexCache;
179
180 typedef typename RowCache::LocalFunctionSpace LFSV;
181 typedef typename ColCache::LocalFunctionSpace LFSU;
182
183 typedef typename LFSV::Traits::DOFIndex RowDOFIndex;
184 typedef typename LFSV::Traits::GridFunctionSpace::Ordering::Traits::ContainerIndex RowContainerIndex;
185
186 typedef typename LFSU::Traits::DOFIndex ColDOFIndex;
187 typedef typename LFSU::Traits::GridFunctionSpace::Ordering::Traits::ContainerIndex ColContainerIndex;
188
189 using BaseT::rowIndexCache;
190 using BaseT::colIndexCache;
191 using BaseT::N;
192 using BaseT::M;
193
194 // Explicitly pull in operator() from the base class to work around a problem
195 // with clang not finding the const overloads of the operator from the base class.
196 using BaseT::operator();
197
198 UncachedMatrixView()
199 {}
200
201 UncachedMatrixView(Container& container)
202 : BaseT(container)
203 {}
204
205 void commit()
206 {}
207
208 template<typename LC>
209 void write(const LC& local_container)
210 {
211 for (size_type i = 0; i < N(); ++i)
212 for (size_type j = 0; j < M(); ++j)
213 container()(rowIndexCache().containerIndex(i),colIndexCache().containerIndex(j)) = local_container.getEntry(i,j);
214 }
215
216 template<typename LC>
217 void add(const LC& local_container)
218 {
219 for (size_type i = 0; i < N(); ++i)
220 for (size_type j = 0; j < M(); ++j)
221 container()(rowIndexCache().containerIndex(i),colIndexCache().containerIndex(j)) += local_container.getEntry(i,j);
222 }
223
224
225
226 ElementType& operator()(size_type i, size_type j)
227 {
228 return container()(rowIndexCache().containerIndex(i),colIndexCache().containerIndex(j));
229 }
230
231 // disable this function if DOFIndex and ContainerIndex have the same type - required for interoperability
232 // with function spaces based on dune-functions bases
233 template<typename RDI, typename CDI>
234 std::enable_if_t<
235 (std::is_same<RDI,RowDOFIndex>{} and std::is_same<CDI,ColDOFIndex>{} and not
236 (std::is_same<RDI,RowContainerIndex>{} and std::is_same<CDI,ColContainerIndex>{})),
237 const ElementType&
238 >
239 operator()(const RDI& i, const CDI& j)
240 {
241 return container()(rowIndexCache().containerIndex(i),colIndexCache().containerIndex(j));
242 }
243
244 ElementType& operator()(const RowContainerIndex& i, const ColContainerIndex& j)
245 {
246 return container()(i,j);
247 }
248
249 ElementType& operator()(const RowContainerIndex& i, size_type j)
250 {
251 return container()(i,colIndexCache().containerIndex(j));
252 }
253
254 ElementType& operator()(size_type i, const ColContainerIndex& j)
255 {
256 return container()(rowIndexCache().containerIndex(i),j);
257 }
258
259 void add(size_type i, size_type j, const ElementType& v)
260 {
261 container()(rowIndexCache().containerIndex(i),colIndexCache().containerIndex(j)) += v;
262 }
263
264 // disable this function if DOFIndex and ContainerIndex have the same type - required for interoperability
265 // with function spaces based on dune-functions bases
266 template<typename RDI, typename CDI>
267 std::enable_if_t<
268 (std::is_same<RDI,RowDOFIndex>{} and std::is_same<CDI,ColDOFIndex>{} and not
269 (std::is_same<RDI,RowContainerIndex>{} and std::is_same<CDI,ColContainerIndex>{}))
270 >
271 add(const RDI& i, const CDI& j, const ElementType& v)
272 {
273 container()(rowIndexCache().containerIndex(i),colIndexCache().containerIndex(j)) += v;
274 }
275
276 void add(const RowContainerIndex& i, const ColContainerIndex& j, const ElementType& v)
277 {
278 container()(i,j) += v;
279 }
280
281 void add(const RowContainerIndex& i, size_type j, const ElementType& v)
282 {
283 container()(i,colIndexCache().containerIndex(j)) += v;
284 }
285
286 void add(size_type i, const ColContainerIndex& j, const ElementType& v)
287 {
288 container()(rowIndexCache().containerIndex(i),j) += v;
289 }
290
291 Container& container()
292 {
293 return *(this->_container);
294 }
295
296 };
297
298
299 } // namespace PDELab
300} // namespace Dune
301
302#endif // DUNE_PDELAB_BACKEND_COMMON_UNCACHEDMATRIXVIEW_HH
Dune namespace.
Definition: alignedallocator.hh:13
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jul 15, 22:36, 2024)