Dune Core Modules (2.9.1)

multitypeblockvector.hh
1// SPDX-FileCopyrightText: Copyright (C) DUNE Project contributors, see file LICENSE.md in module root
2// SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception
3// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
4// vi: set et ts=4 sw=2 sts=2:
5#ifndef DUNE_ISTL_MULTITYPEBLOCKVECTOR_HH
6#define DUNE_ISTL_MULTITYPEBLOCKVECTOR_HH
7
8#include <cmath>
9#include <iostream>
10#include <tuple>
11
14#include <dune/common/hybridutilities.hh>
16#include <dune/common/std/type_traits.hh>
17
18#include "istlexception.hh"
19
20// forward declaration
21namespace Dune {
22 template < typename... Args >
23 class MultiTypeBlockVector;
24}
25
26#include "gsetc.hh"
27
28namespace Dune {
29
41 template <typename... Args>
42 struct FieldTraits< MultiTypeBlockVector<Args...> >
43 {
44 using field_type = typename MultiTypeBlockVector<Args...>::field_type;
45 using real_type = typename FieldTraits<field_type>::real_type;
46 };
56 template < typename... Args >
58 : public std::tuple<Args...>
59 {
61 typedef std::tuple<Args...> TupleType;
62 public:
63
65 using size_type = std::size_t;
66
70 using TupleType::TupleType;
71
75 typedef MultiTypeBlockVector<Args...> type;
76
83
84 // make sure that we have an std::common_type: using an assertion produces a more readable error message
85 // than a compiler template instantiation error
86 static_assert ( sizeof...(Args) == 0 or not std::is_same_v<field_type, Std::nonesuch>,
87 "No std::common_type implemented for the given field_types of the Args. Please provide an implementation and check that a FieldTraits class is present for your type.");
88
89
95 static constexpr size_type size()
96 {
97 return sizeof...(Args);
98 }
99
102 static constexpr size_type N()
103 {
104 return sizeof...(Args);
105 }
106
113 [[deprecated("Use method 'N' instead")]]
114 int count() const
115 {
116 return sizeof...(Args);
117 }
118
121 {
122 size_type result = 0;
123 Hybrid::forEach(std::make_index_sequence<N()>{},
124 [&](auto i){result += std::get<i>(*this).dim();});
125
126 return result;
127 }
128
147 template< size_type index >
148 typename std::tuple_element<index,TupleType>::type&
149 operator[] ([[maybe_unused]] const std::integral_constant< size_type, index > indexVariable)
150 {
151 return std::get<index>(*this);
152 }
153
159 template< size_type index >
160 const typename std::tuple_element<index,TupleType>::type&
161 operator[] ([[maybe_unused]] const std::integral_constant< size_type, index > indexVariable) const
162 {
163 return std::get<index>(*this);
164 }
165
168 template<typename T>
169 void operator= (const T& newval) {
170 Dune::Hybrid::forEach(*this, [&](auto&& entry) {
171 entry = newval;
172 });
173 }
174
178 void operator+= (const type& newv) {
179 using namespace Dune::Hybrid;
180 forEach(integralRange(Hybrid::size(*this)), [&](auto&& i) {
181 (*this)[i] += newv[i];
182 });
183 }
184
188 void operator-= (const type& newv) {
189 using namespace Dune::Hybrid;
190 forEach(integralRange(Hybrid::size(*this)), [&](auto&& i) {
191 (*this)[i] -= newv[i];
192 });
193 }
194
196 template<class T,
197 std::enable_if_t< IsNumber<T>::value, int> = 0>
198 void operator*= (const T& w) {
199 Hybrid::forEach(*this, [&](auto&& entry) {
200 entry *= w;
201 });
202 }
203
205 template<class T,
206 std::enable_if_t< IsNumber<T>::value, int> = 0>
207 void operator/= (const T& w) {
208 Hybrid::forEach(*this, [&](auto&& entry) {
209 entry /= w;
210 });
211 }
212
213 field_type operator* (const type& newv) const {
214 using namespace Dune::Hybrid;
215 return accumulate(integralRange(Hybrid::size(*this)), field_type(0), [&](auto&& a, auto&& i) {
216 return a + (*this)[i]*newv[i];
217 });
218 }
219
220 field_type dot (const type& newv) const {
221 using namespace Dune::Hybrid;
222 return accumulate(integralRange(Hybrid::size(*this)), field_type(0), [&](auto&& a, auto&& i) {
223 return a + (*this)[i].dot(newv[i]);
224 });
225 }
226
229 auto one_norm() const {
230 using namespace Dune::Hybrid;
231 return accumulate(*this, typename FieldTraits<field_type>::real_type(0), [&](auto&& a, auto&& entry) {
232 return a + entry.one_norm();
233 });
234 }
235
238 auto one_norm_real() const {
239 using namespace Dune::Hybrid;
240 return accumulate(*this, typename FieldTraits<field_type>::real_type(0), [&](auto&& a, auto&& entry) {
241 return a + entry.one_norm_real();
242 });
243 }
244
247 typename FieldTraits<field_type>::real_type two_norm2() const {
248 using namespace Dune::Hybrid;
249 return accumulate(*this, typename FieldTraits<field_type>::real_type(0), [&](auto&& a, auto&& entry) {
250 return a + entry.two_norm2();
251 });
252 }
253
256 typename FieldTraits<field_type>::real_type two_norm() const {return sqrt(this->two_norm2());}
257
260 typename FieldTraits<field_type>::real_type infinity_norm() const
261 {
262 using namespace Dune::Hybrid;
263 using std::max;
264 using real_type = typename FieldTraits<field_type>::real_type;
265
266 real_type result = 0.0;
267 // Compute max norm tracking appearing nan values
268 // if the field type supports nan.
269 if constexpr (HasNaN<field_type>()) {
270 // This variable will preserve any nan value
271 real_type nanTracker = 1.0;
272 using namespace Dune::Hybrid; // needed for icc, see issue #31
273 forEach(*this, [&](auto&& entry) {
274 real_type entryNorm = entry.infinity_norm();
275 result = max(entryNorm, result);
276 nanTracker += entryNorm;
277 });
278 // Incorporate possible nan value into result
279 result *= (nanTracker / nanTracker);
280 } else {
281 using namespace Dune::Hybrid; // needed for icc, see issue #31
282 forEach(*this, [&](auto&& entry) {
283 result = max(entry.infinity_norm(), result);
284 });
285 }
286 return result;
287 }
288
291 typename FieldTraits<field_type>::real_type infinity_norm_real() const
292 {
293 using namespace Dune::Hybrid;
294 using std::max;
295 using real_type = typename FieldTraits<field_type>::real_type;
296
297 real_type result = 0.0;
298 // Compute max norm tracking appearing nan values
299 // if the field type supports nan.
300 if constexpr (HasNaN<field_type>()) {
301 // This variable will preserve any nan value
302 real_type nanTracker = 1.0;
303 using namespace Dune::Hybrid; // needed for icc, see issue #31
304 forEach(*this, [&](auto&& entry) {
305 real_type entryNorm = entry.infinity_norm_real();
306 result = max(entryNorm, result);
307 nanTracker += entryNorm;
308 });
309 // Incorporate possible nan value into result
310 result *= (nanTracker / nanTracker);
311 } else {
312 using namespace Dune::Hybrid; // needed for icc, see issue #31
313 forEach(*this, [&](auto&& entry) {
314 result = max(entry.infinity_norm_real(), result);
315 });
316 }
317 return result;
318 }
319
324 template<typename Ta>
325 void axpy (const Ta& a, const type& y) {
326 using namespace Dune::Hybrid;
327 forEach(integralRange(Hybrid::size(*this)), [&](auto&& i) {
328 (*this)[i].axpy(a, y[i]);
329 });
330 }
331
332 };
333
334
335
338 template <typename... Args>
339 std::ostream& operator<< (std::ostream& s, const MultiTypeBlockVector<Args...>& v) {
340 using namespace Dune::Hybrid;
341 forEach(integralRange(Dune::Hybrid::size(v)), [&](auto&& i) {
342 s << "\t(" << i << "):\n" << v[i] << "\n";
343 });
344 return s;
345 }
346
347} // end namespace Dune
348
349namespace std
350{
355 template <size_t i, typename... Args>
356 struct tuple_element<i,Dune::MultiTypeBlockVector<Args...> >
357 {
358 using type = typename std::tuple_element<i, std::tuple<Args...> >::type;
359 };
360}
361
362#endif
A Vector class to support different block types.
Definition: multitypeblockvector.hh:59
Provides the functions dot(a,b) := and dotT(a,b) := .
Type traits to determine the type of reals (when working with complex numbers)
typename detected_or< nonesuch, Op, Args... >::type detected_t
Returns Op<Args...> if that is valid; otherwise returns nonesuch.
Definition: type_traits.hh:170
constexpr void forEach(Range &&range, F &&f)
Range based for loop.
Definition: hybridutilities.hh:268
constexpr auto integralRange(const Begin &begin, const End &end)
Create an integral range.
Definition: hybridutilities.hh:184
constexpr T accumulate(Range &&range, T value, F &&f)
Accumulate values.
Definition: hybridutilities.hh:291
void operator=(const T &newval)
Assignment operator.
Definition: multitypeblockvector.hh:169
std::size_t size_type
Type used for vector sizes.
Definition: multitypeblockvector.hh:65
int count() const
Definition: multitypeblockvector.hh:114
static constexpr size_type N()
Number of elements.
Definition: multitypeblockvector.hh:102
Std::detected_t< std::common_type_t, typename FieldTraits< std::decay_t< Args > >::field_type... > field_type
The type used for scalars.
Definition: multitypeblockvector.hh:82
static constexpr size_type size()
Return the number of non-zero vector entries.
Definition: multitypeblockvector.hh:95
std::tuple_element< index, TupleType >::type & operator[](const std::integral_constant< size_type, index > indexVariable)
Random-access operator.
Definition: multitypeblockvector.hh:149
FieldTraits< field_type >::real_type two_norm() const
Compute the Euclidean norm.
Definition: multitypeblockvector.hh:256
size_type dim() const
Number of scalar elements.
Definition: multitypeblockvector.hh:120
void operator*=(const T &w)
Multiplication with a scalar.
Definition: multitypeblockvector.hh:198
void axpy(const Ta &a, const type &y)
Axpy operation on this vector (*this += a * y)
Definition: multitypeblockvector.hh:325
void operator/=(const T &w)
Division by a scalar.
Definition: multitypeblockvector.hh:207
MultiTypeBlockVector< Args... > type
Definition: multitypeblockvector.hh:75
FieldTraits< field_type >::real_type infinity_norm_real() const
Compute the simplified maximum norm (uses 1-norm for complex values)
Definition: multitypeblockvector.hh:291
auto one_norm() const
Compute the 1-norm.
Definition: multitypeblockvector.hh:229
void operator-=(const type &newv)
Definition: multitypeblockvector.hh:188
void operator+=(const type &newv)
Definition: multitypeblockvector.hh:178
FieldTraits< field_type >::real_type infinity_norm() const
Compute the maximum norm.
Definition: multitypeblockvector.hh:260
auto one_norm_real() const
Compute the simplified 1-norm (uses 1-norm also for complex values)
Definition: multitypeblockvector.hh:238
FieldTraits< field_type >::real_type two_norm2() const
Compute the squared Euclidean norm.
Definition: multitypeblockvector.hh:247
auto max(ADLTag< 0 >, const V &v1, const V &v2)
implements binary Simd::max()
Definition: defaults.hh:81
Simple iterative methods like Jacobi, Gauss-Seidel, SOR, SSOR, etc. in a generic way.
Dune namespace.
Definition: alignedallocator.hh:13
STL namespace.
Whether this type has a value of NaN.
Definition: typetraits.hh:212
Traits for type conversions and type information.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jul 15, 22:36, 2024)