Dune Core Modules (unstable)

multitypeblockvector.hh
1 // SPDX-FileCopyrightText: Copyright © 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 
13 #include <dune/common/ftraits.hh>
14 #include <dune/common/hybridutilities.hh>
16 #include <dune/common/std/type_traits.hh>
17 
18 #include "istlexception.hh"
19 
20 // forward declaration
21 namespace Dune {
22  template < typename... Args >
23  class MultiTypeBlockVector;
24 }
25 
26 #include "gsetc.hh"
27 
28 namespace 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 
108  size_type dim() const
109  {
110  size_type result = 0;
111  Hybrid::forEach(std::make_index_sequence<N()>{},
112  [&](auto i){result += std::get<i>(*this).dim();});
113 
114  return result;
115  }
116 
135  template< size_type index >
136  typename std::tuple_element<index,TupleType>::type&
137  operator[] ([[maybe_unused]] const std::integral_constant< size_type, index > indexVariable)
138  {
139  return std::get<index>(*this);
140  }
141 
147  template< size_type index >
148  const typename std::tuple_element<index,TupleType>::type&
149  operator[] ([[maybe_unused]] const std::integral_constant< size_type, index > indexVariable) const
150  {
151  return std::get<index>(*this);
152  }
153 
156  template<typename T>
157  void operator= (const T& newval) {
158  Dune::Hybrid::forEach(*this, [&](auto&& entry) {
159  entry = newval;
160  });
161  }
162 
166  void operator+= (const type& newv) {
167  using namespace Dune::Hybrid;
168  forEach(integralRange(Hybrid::size(*this)), [&](auto&& i) {
169  (*this)[i] += newv[i];
170  });
171  }
172 
176  void operator-= (const type& newv) {
177  using namespace Dune::Hybrid;
178  forEach(integralRange(Hybrid::size(*this)), [&](auto&& i) {
179  (*this)[i] -= newv[i];
180  });
181  }
182 
184  template<class T,
185  std::enable_if_t< IsNumber<T>::value, int> = 0>
186  void operator*= (const T& w) {
187  Hybrid::forEach(*this, [&](auto&& entry) {
188  entry *= w;
189  });
190  }
191 
193  template<class T,
194  std::enable_if_t< IsNumber<T>::value, int> = 0>
195  void operator/= (const T& w) {
196  Hybrid::forEach(*this, [&](auto&& entry) {
197  entry /= w;
198  });
199  }
200 
201  field_type operator* (const type& newv) const {
202  using namespace Dune::Hybrid;
203  return accumulate(integralRange(Hybrid::size(*this)), field_type(0), [&](auto&& a, auto&& i) {
204  return a + (*this)[i]*newv[i];
205  });
206  }
207 
208  field_type dot (const type& newv) const {
209  using namespace Dune::Hybrid;
210  return accumulate(integralRange(Hybrid::size(*this)), field_type(0), [&](auto&& a, auto&& i) {
211  return a + (*this)[i].dot(newv[i]);
212  });
213  }
214 
217  auto one_norm() const {
218  using namespace Dune::Hybrid;
219  return accumulate(*this, typename FieldTraits<field_type>::real_type(0), [&](auto&& a, auto&& entry) {
220  return a + entry.one_norm();
221  });
222  }
223 
226  auto one_norm_real() const {
227  using namespace Dune::Hybrid;
228  return accumulate(*this, typename FieldTraits<field_type>::real_type(0), [&](auto&& a, auto&& entry) {
229  return a + entry.one_norm_real();
230  });
231  }
232 
235  typename FieldTraits<field_type>::real_type two_norm2() const {
236  using namespace Dune::Hybrid;
237  return accumulate(*this, typename FieldTraits<field_type>::real_type(0), [&](auto&& a, auto&& entry) {
238  return a + entry.two_norm2();
239  });
240  }
241 
244  typename FieldTraits<field_type>::real_type two_norm() const {return sqrt(this->two_norm2());}
245 
248  typename FieldTraits<field_type>::real_type infinity_norm() const
249  {
250  using namespace Dune::Hybrid;
251  using std::max;
252  using real_type = typename FieldTraits<field_type>::real_type;
253 
254  real_type result = 0.0;
255  // Compute max norm tracking appearing nan values
256  // if the field type supports nan.
257  if constexpr (HasNaN<field_type>()) {
258  // This variable will preserve any nan value
259  real_type nanTracker = 1.0;
260  using namespace Dune::Hybrid; // needed for icc, see issue #31
261  forEach(*this, [&](auto&& entry) {
262  real_type entryNorm = entry.infinity_norm();
263  result = max(entryNorm, result);
264  nanTracker += entryNorm;
265  });
266  // Incorporate possible nan value into result
267  result *= (nanTracker / nanTracker);
268  } else {
269  using namespace Dune::Hybrid; // needed for icc, see issue #31
270  forEach(*this, [&](auto&& entry) {
271  result = max(entry.infinity_norm(), result);
272  });
273  }
274  return result;
275  }
276 
279  typename FieldTraits<field_type>::real_type infinity_norm_real() const
280  {
281  using namespace Dune::Hybrid;
282  using std::max;
283  using real_type = typename FieldTraits<field_type>::real_type;
284 
285  real_type result = 0.0;
286  // Compute max norm tracking appearing nan values
287  // if the field type supports nan.
288  if constexpr (HasNaN<field_type>()) {
289  // This variable will preserve any nan value
290  real_type nanTracker = 1.0;
291  using namespace Dune::Hybrid; // needed for icc, see issue #31
292  forEach(*this, [&](auto&& entry) {
293  real_type entryNorm = entry.infinity_norm_real();
294  result = max(entryNorm, result);
295  nanTracker += entryNorm;
296  });
297  // Incorporate possible nan value into result
298  result *= (nanTracker / nanTracker);
299  } else {
300  using namespace Dune::Hybrid; // needed for icc, see issue #31
301  forEach(*this, [&](auto&& entry) {
302  result = max(entry.infinity_norm_real(), result);
303  });
304  }
305  return result;
306  }
307 
312  template<typename Ta>
313  void axpy (const Ta& a, const type& y) {
314  using namespace Dune::Hybrid;
315  forEach(integralRange(Hybrid::size(*this)), [&](auto&& i) {
316  (*this)[i].axpy(a, y[i]);
317  });
318  }
319 
320  };
321 
322 
323 
326  template <typename... Args>
327  std::ostream& operator<< (std::ostream& s, const MultiTypeBlockVector<Args...>& v) {
328  using namespace Dune::Hybrid;
329  forEach(integralRange(Dune::Hybrid::size(v)), [&](auto&& i) {
330  s << "\t(" << i << "):\n" << v[i] << "\n";
331  });
332  return s;
333  }
334 
335 } // end namespace Dune
336 
337 namespace std
338 {
343  template <size_t i, typename... Args>
344  struct tuple_element<i,Dune::MultiTypeBlockVector<Args...> >
345  {
346  using type = typename std::tuple_element<i, std::tuple<Args...> >::type;
347  };
348 
353  template <typename... Args>
354  struct tuple_size<Dune::MultiTypeBlockVector<Args...> >
355  : std::integral_constant<std::size_t, sizeof...(Args)>
356  {};
357 }
358 
359 #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:173
constexpr auto size(const T &t)
Size query.
Definition: hybridutilities.hh:73
constexpr void forEach(Range &&range, F &&f)
Range based for loop.
Definition: hybridutilities.hh:256
constexpr auto max
Function object that returns the greater of the given values.
Definition: hybridutilities.hh:484
constexpr auto integralRange(const Begin &begin, const End &end)
Create an integral range.
Definition: hybridutilities.hh:172
constexpr T accumulate(Range &&range, T value, F &&f)
Accumulate values.
Definition: hybridutilities.hh:279
void operator=(const T &newval)
Assignment operator.
Definition: multitypeblockvector.hh:157
std::size_t size_type
Type used for vector sizes.
Definition: multitypeblockvector.hh:65
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[]([[maybe_unused]] const std::integral_constant< size_type, index > indexVariable)
Random-access operator.
Definition: multitypeblockvector.hh:137
size_type dim() const
Number of scalar elements.
Definition: multitypeblockvector.hh:108
void operator*=(const T &w)
Multiplication with a scalar.
Definition: multitypeblockvector.hh:186
void axpy(const Ta &a, const type &y)
Axpy operation on this vector (*this += a * y)
Definition: multitypeblockvector.hh:313
void operator/=(const T &w)
Division by a scalar.
Definition: multitypeblockvector.hh:195
MultiTypeBlockVector< Args... > type
Definition: multitypeblockvector.hh:75
auto one_norm() const
Compute the 1-norm.
Definition: multitypeblockvector.hh:217
void operator-=(const type &newv)
Definition: multitypeblockvector.hh:176
FieldTraits< field_type >::real_type infinity_norm_real() const
Compute the simplified maximum norm (uses 1-norm for complex values)
Definition: multitypeblockvector.hh:279
FieldTraits< field_type >::real_type infinity_norm() const
Compute the maximum norm.
Definition: multitypeblockvector.hh:248
void operator+=(const type &newv)
Definition: multitypeblockvector.hh:166
FieldTraits< field_type >::real_type two_norm2() const
Compute the squared Euclidean norm.
Definition: multitypeblockvector.hh:235
FieldTraits< field_type >::real_type two_norm() const
Compute the Euclidean norm.
Definition: multitypeblockvector.hh:244
auto one_norm_real() const
Compute the simplified 1-norm (uses 1-norm also for complex values)
Definition: multitypeblockvector.hh:226
Simple iterative methods like Jacobi, Gauss-Seidel, SOR, SSOR, etc. in a generic way.
Dune namespace.
Definition: alignedallocator.hh:13
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.80.0 (Apr 27, 22:29, 2024)