DUNE PDELab (2.7)

vector.hh
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_PDELAB_BACKEND_SIMPLE_VECTOR_HH
4#define DUNE_PDELAB_BACKEND_SIMPLE_VECTOR_HH
5
6#include <algorithm>
7#include <functional>
8#include <memory>
9#include <numeric>
10
13#include <dune/istl/bvector.hh>
14
15#include <dune/pdelab/gridfunctionspace/gridfunctionspace.hh>
16#include <dune/pdelab/gridfunctionspace/lfsindexcache.hh>
18#include <dune/pdelab/backend/interface.hh>
19#include <dune/pdelab/backend/common/uncachedvectorview.hh>
20#include <dune/pdelab/backend/simple/descriptors.hh>
21
22namespace Dune {
23 namespace PDELab {
24
25 namespace Simple {
26
27 namespace {
28
29 // For some reason std::bind cannot directly deduce the correct version
30 // of Dune::fvmeta::abs2, so we package it into a functor to help it along.
31 template<typename K>
32 struct abs2
33 {
34 typename FieldTraits<K>::real_type operator()(const K& k) const
35 {
36 return Dune::fvmeta::abs2(k);
37 }
38 };
39
40 }
41
42 template<typename GFS, typename C>
43 class VectorContainer
44 : public Backend::impl::Wrapper<C>
45 {
46
47 friend Backend::impl::Wrapper<C>;
48
49 public:
50 typedef C Container;
51 typedef typename Container::value_type ElementType;
52 typedef ElementType E;
53
54 // for ISTL solver compatibility
55 typedef ElementType field_type;
56
57 typedef GFS GridFunctionSpace;
58 typedef typename Container::size_type size_type;
59
60 typedef typename GFS::Ordering::Traits::ContainerIndex ContainerIndex;
61
62 typedef typename Container::iterator iterator;
63 typedef typename Container::const_iterator const_iterator;
64
65 template<typename LFSCache>
66 using LocalView = UncachedVectorView<VectorContainer,LFSCache>;
67
68 template<typename LFSCache>
69 using ConstLocalView = ConstUncachedVectorView<const VectorContainer,LFSCache>;
70
71
72 VectorContainer(const VectorContainer& rhs)
73 : _gfs(rhs._gfs)
74 , _container(std::make_shared<Container>(*(rhs._container)))
75 {}
76
77 VectorContainer (const GFS& gfs, Backend::attached_container = Backend::attached_container())
78 : _gfs(gfs)
79 , _container(std::make_shared<Container>(gfs.ordering().blockCount()))
80 {}
81
83 VectorContainer(const GFS& gfs, Backend::unattached_container)
84 : _gfs(gfs)
85 {}
86
92 VectorContainer (const GFS& gfs, Container& container)
93 : _gfs(gfs)
94 , _container(stackobject_to_shared_ptr(container))
95 {
96 _container->resize(gfs.ordering().blockCount());
97 }
98
99 VectorContainer (const GFS& gfs, const E& e)
100 : _gfs(gfs)
101 , _container(std::make_shared<Container>(gfs.ordering().blockCount(),e))
102 {}
103
104 void detach()
105 {
106 _container.reset();
107 }
108
109 void attach(std::shared_ptr<Container> container)
110 {
111 _container = container;
112 }
113
114 bool attached() const
115 {
116 return bool(_container);
117 }
118
119 const std::shared_ptr<Container>& storage() const
120 {
121 return _container;
122 }
123
124 size_type N() const
125 {
126 return _container->size();
127 }
128
129 VectorContainer& operator=(const VectorContainer& r)
130 {
131 if (this == &r)
132 return *this;
133 if (attached())
134 {
135 (*_container) = (*r._container);
136 }
137 else
138 {
139 _container = std::make_shared<Container>(*(r._container));
140 }
141 return *this;
142 }
143
144 VectorContainer& operator=(const E& e)
145 {
146 std::fill(_container->begin(),_container->end(),e);
147 return *this;
148 }
149
150 VectorContainer& operator*=(const E& e)
151 {
152 std::transform(_container->begin(),_container->end(),_container->begin(),
153 std::bind(std::multiplies<E>(),e,std::placeholders::_1));
154 return *this;
155 }
156
157
158 VectorContainer& operator+=(const E& e)
159 {
160 std::transform(_container->begin(),_container->end(),_container->begin(),
161 std::bind(std::plus<E>(),e,std::placeholders::_1));
162 return *this;
163 }
164
165 VectorContainer& operator+=(const VectorContainer& y)
166 {
167 std::transform(_container->begin(),_container->end(),y._container->begin(),
168 _container->begin(),std::plus<E>());
169 return *this;
170 }
171
172 VectorContainer& operator-= (const VectorContainer& y)
173 {
174 std::transform(_container->begin(),_container->end(),y._container->begin(),
175 _container->begin(),std::minus<E>());
176 return *this;
177 }
178
179 E& operator[](const ContainerIndex& ci)
180 {
181 return (*_container)[ci[0]];
182 }
183
184 const E& operator[](const ContainerIndex& ci) const
185 {
186 return (*_container)[ci[0]];
187 }
188
189 typename Dune::template FieldTraits<E>::real_type two_norm() const
190 {
191 using namespace std::placeholders;
192 typedef typename Dune::template FieldTraits<E>::real_type Real;
193 return std::sqrt(std::accumulate(_container->begin(),_container->end(),Real(0),std::bind(std::plus<Real>(),_1,std::bind(abs2<E>(),_2))));
194 }
195
196 typename Dune::template FieldTraits<E>::real_type one_norm() const
197 {
198 typedef typename Dune::template FieldTraits<E>::real_type Real;
199 return std::accumulate(_container->begin(),_container->end(),Real(0),
200 [](const auto& n, const auto& e) {
201 using std::abs;
202 return n + abs(e);
203 });
204 }
205
206 typename Dune::template FieldTraits<E>::real_type infinity_norm() const
207 {
208 if (_container->size() == 0)
209 return 0;
210 using std::abs;
211 return abs(*std::max_element(_container->begin(),_container->end(),
212 [](const auto& a, const auto& b) {
213 using std::abs;
214 return abs(a) < abs(b);
215 }));
216 }
217
218 E operator*(const VectorContainer& y) const
219 {
220 return std::inner_product(_container->begin(),_container->end(),y._container->begin(),E(0));
221 }
222
223 E dot(const VectorContainer& y) const
224 {
225 return std::inner_product(_container->begin(),_container->end(),y._container->begin(),E(0),std::plus<E>(),Dune::dot<E,E>);
226 }
227
228 VectorContainer& axpy(const E& a, const VectorContainer& y)
229 {
230 using namespace std::placeholders;
231 std::transform(_container->begin(),_container->end(),y._container->begin(),
232 _container->begin(),std::bind(std::plus<E>(),_1,std::bind(std::multiplies<E>(),a,_2)));
233 return *this;
234 }
235
236 // for debugging and AMG access
237 Container& base ()
238 {
239 return *_container;
240 }
241
242 const Container& base () const
243 {
244 return *_container;
245 }
246
247 private:
248
249 Container& native ()
250 {
251 return *_container;
252 }
253
254 const Container& native () const
255 {
256 return *_container;
257 }
258
259 public:
260
261 iterator begin()
262 {
263 return _container->begin();
264 }
265
266 const_iterator begin() const
267 {
268 return _container->begin();
269 }
270
271 iterator end()
272 {
273 return _container->end();
274 }
275
276 const_iterator end() const
277 {
278 return _container->end();
279 }
280
281 size_t flatsize() const
282 {
283 return _container->size();
284 }
285
286 const GFS& gridFunctionSpace() const
287 {
288 return _gfs;
289 }
290
291 private:
292 const GFS& _gfs;
293 std::shared_ptr<Container> _container;
294 };
295
296
297 } // namespace Simple
298
299
300#ifndef DOXYGEN
301
302 template<typename GFS, typename E>
303 struct SimpleVectorSelectorHelper
304 {
305
306 using vector_type = typename GFS::Traits::Backend::template vector_type<E>;
307
308 using Type = Simple::VectorContainer<GFS,vector_type>;
309
310 };
311
312 namespace Backend {
313 namespace impl {
314
315 template<template<typename> class Container, typename GFS, typename E>
316 struct BackendVectorSelectorHelper<Simple::VectorBackend<Container>, GFS, E>
317 : public SimpleVectorSelectorHelper<GFS,E>
318 {};
319
320 } // namespace impl
321 } // namespace Backend
322
323#endif // DOXYGEN
324
325 } // namespace PDELab
326} // namespace Dune
327
328#endif // DUNE_PDELAB_BACKEND_SIMPLE_VECTOR_HH
Various tags for influencing backend behavior.
This file implements a vector space as a tensor product of a given vector space. The number of compon...
Implements a vector constructed from a given type representing a field and a compile-time given size.
auto dot(const A &a, const B &b) -> typename std::enable_if<!IsVector< A >::value &&!std::is_same< typename FieldTraits< A >::field_type, typename FieldTraits< A >::real_type > ::value, decltype(conj(a) *b)>::type
computes the dot product for fundamental data types according to Petsc's VectDot function: dot(a,...
Definition: dotproduct.hh:40
constexpr index_constant< 1 > _1
Compile time index with value 1.
Definition: indices.hh:54
T accumulate(Range &&range, T value, F &&f)
Accumulate values.
Definition: hybridutilities.hh:290
Dune namespace.
Definition: alignedallocator.hh:14
shared_ptr< T > stackobject_to_shared_ptr(T &t)
Create a shared_ptr for a stack-allocated object.
Definition: shared_ptr.hh:75
STL namespace.
This file implements the class shared_ptr (a reference counting pointer), for those systems that don'...
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jul 15, 22:36, 2024)