Dune Core Modules (unstable)

multitypeblockmatrix.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_MULTITYPEBLOCKMATRIX_HH
6 #define DUNE_ISTL_MULTITYPEBLOCKMATRIX_HH
7 
8 #include <cmath>
9 #include <iostream>
10 #include <tuple>
11 
12 #include <dune/common/hybridutilities.hh>
13 
14 #include "istlexception.hh"
15 
16 // forward declaration
17 namespace Dune
18 {
19  template<typename FirstRow, typename... Args>
20  class MultiTypeBlockMatrix;
21 
22  template<int I, int crow, int remain_row>
23  class MultiTypeBlockMatrix_Solver;
24 }
25 
26 #include "gsetc.hh"
27 
28 namespace Dune {
29 
43  template<typename FirstRow, typename... Args>
45  : public std::tuple<FirstRow, Args...>
46  {
47  using ParentType = std::tuple<FirstRow, Args...>;
48  public:
49 
51  using ParentType::ParentType;
52 
56  typedef MultiTypeBlockMatrix<FirstRow, Args...> type;
57 
59  using size_type = std::size_t;
60 
61  typedef typename FirstRow::field_type field_type;
62 
64  static constexpr size_type N()
65  {
66  return 1+sizeof...(Args);
67  }
68 
70  static constexpr size_type M()
71  {
72  return FirstRow::size();
73  }
74 
91  template< size_type index >
92  auto
93  operator[] ([[maybe_unused]] const std::integral_constant< size_type, index > indexVariable)
94  -> decltype(std::get<index>(*this))
95  {
96  return std::get<index>(*this);
97  }
98 
104  template< size_type index >
105  auto
106  operator[] ([[maybe_unused]] const std::integral_constant< size_type, index > indexVariable) const
107  -> decltype(std::get<index>(*this))
108  {
109  return std::get<index>(*this);
110  }
111 
115  template<typename T>
116  void operator= (const T& newval) {
117  using namespace Dune::Hybrid;
118  auto size = index_constant<1+sizeof...(Args)>();
119  // Since Dune::Hybrid::size(MultiTypeBlockMatrix) is not implemented,
120  // we cannot use a plain forEach(*this, ...). This could be achieved,
121  // e.g., by implementing a static size() method.
122  forEach(integralRange(size), [&](auto&& i) {
123  (*this)[i] = newval;
124  });
125  }
126 
127  //===== vector space arithmetic
128 
130  MultiTypeBlockMatrix& operator*= (const field_type& k)
131  {
132  auto size = index_constant<N()>();
133  Hybrid::forEach(Hybrid::integralRange(size), [&](auto&& i) {
134  (*this)[i] *= k;
135  });
136 
137  return *this;
138  }
139 
141  MultiTypeBlockMatrix& operator/= (const field_type& k)
142  {
143  auto size = index_constant<N()>();
144  Hybrid::forEach(Hybrid::integralRange(size), [&](auto&& i) {
145  (*this)[i] /= k;
146  });
147 
148  return *this;
149  }
150 
151 
158  {
159  auto size = index_constant<N()>();
160  Hybrid::forEach(Hybrid::integralRange(size), [&](auto&& i) {
161  (*this)[i] += b[i];
162  });
163 
164  return *this;
165  }
166 
173  {
174  auto size = index_constant<N()>();
175  Hybrid::forEach(Hybrid::integralRange(size), [&](auto&& i) {
176  (*this)[i] -= b[i];
177  });
178 
179  return *this;
180  }
181 
184  template<typename X, typename Y>
185  void mv (const X& x, Y& y) const {
186  static_assert(X::size() == M(), "length of x does not match row length");
187  static_assert(Y::size() == N(), "length of y does not match row count");
188  y = 0; //reset y (for mv uses umv)
189  umv(x,y);
190  }
191 
194  template<typename X, typename Y>
195  void umv (const X& x, Y& y) const {
196  static_assert(X::size() == M(), "length of x does not match row length");
197  static_assert(Y::size() == N(), "length of y does not match row count");
198  using namespace Dune::Hybrid;
199  forEach(integralRange(Hybrid::size(y)), [&](auto&& i) {
200  using namespace Dune::Hybrid; // needed for icc, see issue #31
201  forEach(integralRange(Hybrid::size(x)), [&](auto&& j) {
202  (*this)[i][j].umv(x[j], y[i]);
203  });
204  });
205  }
206 
209  template<typename X, typename Y>
210  void mmv (const X& x, Y& y) const {
211  static_assert(X::size() == M(), "length of x does not match row length");
212  static_assert(Y::size() == N(), "length of y does not match row count");
213  using namespace Dune::Hybrid;
214  forEach(integralRange(Hybrid::size(y)), [&](auto&& i) {
215  using namespace Dune::Hybrid; // needed for icc, see issue #31
216  forEach(integralRange(Hybrid::size(x)), [&](auto&& j) {
217  (*this)[i][j].mmv(x[j], y[i]);
218  });
219  });
220  }
221 
224  template<typename AlphaType, typename X, typename Y>
225  void usmv (const AlphaType& alpha, const X& x, Y& y) const {
226  static_assert(X::size() == M(), "length of x does not match row length");
227  static_assert(Y::size() == N(), "length of y does not match row count");
228  using namespace Dune::Hybrid;
229  forEach(integralRange(Hybrid::size(y)), [&](auto&& i) {
230  using namespace Dune::Hybrid; // needed for icc, see issue #31
231  forEach(integralRange(Hybrid::size(x)), [&](auto&& j) {
232  (*this)[i][j].usmv(alpha, x[j], y[i]);
233  });
234  });
235  }
236 
239  template<typename X, typename Y>
240  void mtv (const X& x, Y& y) const {
241  static_assert(X::size() == N(), "length of x does not match number of rows");
242  static_assert(Y::size() == M(), "length of y does not match number of columns");
243  y = 0;
244  umtv(x,y);
245  }
246 
249  template<typename X, typename Y>
250  void umtv (const X& x, Y& y) const {
251  static_assert(X::size() == N(), "length of x does not match number of rows");
252  static_assert(Y::size() == M(), "length of y does not match number of columns");
253  using namespace Dune::Hybrid;
254  forEach(integralRange(Hybrid::size(y)), [&](auto&& i) {
255  using namespace Dune::Hybrid; // needed for icc, see issue #31
256  forEach(integralRange(Hybrid::size(x)), [&](auto&& j) {
257  (*this)[j][i].umtv(x[j], y[i]);
258  });
259  });
260  }
261 
264  template<typename X, typename Y>
265  void mmtv (const X& x, Y& y) const {
266  static_assert(X::size() == N(), "length of x does not match number of rows");
267  static_assert(Y::size() == M(), "length of y does not match number of columns");
268  using namespace Dune::Hybrid;
269  forEach(integralRange(Hybrid::size(y)), [&](auto&& i) {
270  using namespace Dune::Hybrid; // needed for icc, see issue #31
271  forEach(integralRange(Hybrid::size(x)), [&](auto&& j) {
272  (*this)[j][i].mmtv(x[j], y[i]);
273  });
274  });
275  }
276 
279  template<typename X, typename Y>
280  void usmtv (const field_type& alpha, const X& x, Y& y) const {
281  static_assert(X::size() == N(), "length of x does not match number of rows");
282  static_assert(Y::size() == M(), "length of y does not match number of columns");
283  using namespace Dune::Hybrid;
284  forEach(integralRange(Hybrid::size(y)), [&](auto&& i) {
285  using namespace Dune::Hybrid; // needed for icc, see issue #31
286  forEach(integralRange(Hybrid::size(x)), [&](auto&& j) {
287  (*this)[j][i].usmtv(alpha, x[j], y[i]);
288  });
289  });
290  }
291 
294  template<typename X, typename Y>
295  void umhv (const X& x, Y& y) const {
296  static_assert(X::size() == N(), "length of x does not match number of rows");
297  static_assert(Y::size() == M(), "length of y does not match number of columns");
298  using namespace Dune::Hybrid;
299  forEach(integralRange(Hybrid::size(y)), [&](auto&& i) {
300  using namespace Dune::Hybrid; // needed for icc, see issue #31
301  forEach(integralRange(Hybrid::size(x)), [&](auto&& j) {
302  (*this)[j][i].umhv(x[j], y[i]);
303  });
304  });
305  }
306 
309  template<typename X, typename Y>
310  void mmhv (const X& x, Y& y) const {
311  static_assert(X::size() == N(), "length of x does not match number of rows");
312  static_assert(Y::size() == M(), "length of y does not match number of columns");
313  using namespace Dune::Hybrid;
314  forEach(integralRange(Hybrid::size(y)), [&](auto&& i) {
315  using namespace Dune::Hybrid; // needed for icc, see issue #31
316  forEach(integralRange(Hybrid::size(x)), [&](auto&& j) {
317  (*this)[j][i].mmhv(x[j], y[i]);
318  });
319  });
320  }
321 
324  template<typename X, typename Y>
325  void usmhv (const field_type& alpha, const X& x, Y& y) const {
326  static_assert(X::size() == N(), "length of x does not match number of rows");
327  static_assert(Y::size() == M(), "length of y does not match number of columns");
328  using namespace Dune::Hybrid;
329  forEach(integralRange(Hybrid::size(y)), [&](auto&& i) {
330  using namespace Dune::Hybrid; // needed for icc, see issue #31
331  forEach(integralRange(Hybrid::size(x)), [&](auto&& j) {
332  (*this)[j][i].usmhv(alpha, x[j], y[i]);
333  });
334  });
335  }
336 
337 
338  //===== norms
339 
341  auto frobenius_norm2 () const
342  {
343  using field_type_00 = typename std::decay_t<decltype((*this)[Indices::_0][Indices::_0])>::field_type;
344  typename FieldTraits<field_type_00>::real_type sum=0;
345 
346  auto rows = index_constant<N()>();
347  Hybrid::forEach(Hybrid::integralRange(rows), [&](auto&& i) {
349  Hybrid::forEach(Hybrid::integralRange(cols), [&](auto&& j) {
350  sum += (*this)[i][j].frobenius_norm2();
351  });
352  });
353 
354  return sum;
355  }
356 
358  typename FieldTraits<field_type>::real_type frobenius_norm () const
359  {
360  return sqrt(frobenius_norm2());
361  }
362 
364  auto infinity_norm () const
365  {
366  using field_type_00 = typename std::decay_t<decltype((*this)[Indices::_0][Indices::_0])>::field_type;
367  using std::max;
368  typename FieldTraits<field_type_00>::real_type norm=0;
369 
370  auto rows = index_constant<N()>();
371  Hybrid::forEach(Hybrid::integralRange(rows), [&](auto&& i) {
372 
373  typename FieldTraits<field_type_00>::real_type sum=0;
375  Hybrid::forEach(Hybrid::integralRange(cols), [&](auto&& j) {
376  sum += (*this)[i][j].infinity_norm();
377  });
378  norm = max(sum, norm);
379  });
380 
381  return norm;
382  }
383 
385  auto infinity_norm_real () const
386  {
387  using field_type_00 = typename std::decay_t<decltype((*this)[Indices::_0][Indices::_0])>::field_type;
388  using std::max;
389  typename FieldTraits<field_type_00>::real_type norm=0;
390 
391  auto rows = index_constant<N()>();
392  Hybrid::forEach(Hybrid::integralRange(rows), [&](auto&& i) {
393 
394  typename FieldTraits<field_type_00>::real_type sum=0;
396  Hybrid::forEach(Hybrid::integralRange(cols), [&](auto&& j) {
397  sum += (*this)[i][j].infinity_norm_real();
398  });
399  norm = max(sum, norm);
400  });
401 
402  return norm;
403  }
404 
405  };
406 
407 
411  template <typename... Rows>
412  struct FieldTraits< MultiTypeBlockMatrix<Rows...> >
413  {
414  using field_type = typename MultiTypeBlockMatrix<Rows...>::field_type;
415  using real_type = typename FieldTraits<field_type>::real_type;
416  };
417 
418 
424  template<typename T1, typename... Args>
425  std::ostream& operator<< (std::ostream& s, const MultiTypeBlockMatrix<T1,Args...>& m) {
428  using namespace Dune::Hybrid;
429  forEach(integralRange(N), [&](auto&& i) {
430  using namespace Dune::Hybrid; // needed for icc, see issue #31
431  forEach(integralRange(M), [&](auto&& j) {
432  s << "\t(" << i << ", " << j << "): \n" << m[i][j];
433  });
434  });
435  s << std::endl;
436  return s;
437  }
438 
439  //make algmeta_itsteps known
440  template<int I, typename M>
441  struct algmeta_itsteps;
442 
449  template<int I, int crow, int ccol, int remain_col> //MultiTypeBlockMatrix_Solver_Col: iterating over one row
450  class MultiTypeBlockMatrix_Solver_Col { //calculating b- A[i][j]*x[j]
451  public:
455  template <typename Trhs, typename TVector, typename TMatrix, typename K>
456  static void calc_rhs(const TMatrix& A, TVector& x, TVector& v, Trhs& b, const K& w) {
457  std::get<ccol>( std::get<crow>(A) ).mmv( std::get<ccol>(x), b );
459  }
460 
461  };
462  template<int I, int crow, int ccol> //MultiTypeBlockMatrix_Solver_Col recursion end
463  class MultiTypeBlockMatrix_Solver_Col<I,crow,ccol,0> {
464  public:
465  template <typename Trhs, typename TVector, typename TMatrix, typename K>
466  static void calc_rhs(const TMatrix&, TVector&, TVector&, Trhs&, const K&) {}
467  };
468 
469 
470 
477  template<int I, int crow, int remain_row>
479  public:
480 
484  template <typename TVector, typename TMatrix, typename K>
485  static void dbgs(const TMatrix& A, TVector& x, const TVector& b, const K& w) {
486  TVector xold(x);
487  xold=x; //store old x values
489  x *= w;
490  x.axpy(1-w,xold); //improve x
491  }
492  template <typename TVector, typename TMatrix, typename K>
493  static void dbgs(const TMatrix& A, TVector& x, TVector& v, const TVector& b, const K& w) {
494  auto rhs = std::get<crow> (b);
495 
496  MultiTypeBlockMatrix_Solver_Col<I,crow,0, TMatrix::M()>::calc_rhs(A,x,v,rhs,w); // calculate right side of equation
497  //solve on blocklevel I-1
498  using M =
499  typename std::remove_cv<
500  typename std::remove_reference<
501  decltype(std::get<crow>( std::get<crow>(A)))
502  >::type
503  >::type;
504  algmeta_itsteps<I-1,M>::dbgs(std::get<crow>( std::get<crow>(A)), std::get<crow>(x),rhs,w);
506  }
507 
508 
509 
513  template <typename TVector, typename TMatrix, typename K>
514  static void bsorf(const TMatrix& A, TVector& x, const TVector& b, const K& w) {
515  TVector v;
516  v=x; //use latest x values in right side calculation
518 
519  }
520  template <typename TVector, typename TMatrix, typename K> //recursion over all matrix rows (A)
521  static void bsorf(const TMatrix& A, TVector& x, TVector& v, const TVector& b, const K& w) {
522  auto rhs = std::get<crow> (b);
523 
524  MultiTypeBlockMatrix_Solver_Col<I,crow,0,TMatrix::M()>::calc_rhs(A,x,v,rhs,w); // calculate right side of equation
525  //solve on blocklevel I-1
526  using M =
527  typename std::remove_cv<
528  typename std::remove_reference<
529  decltype(std::get<crow>( std::get<crow>(A)))
530  >::type
531  >::type;
532  algmeta_itsteps<I-1,M>::bsorf(std::get<crow>( std::get<crow>(A)), std::get<crow>(v),rhs,w);
533  std::get<crow>(x).axpy(w,std::get<crow>(v));
535  }
536 
540  template <typename TVector, typename TMatrix, typename K>
541  static void bsorb(const TMatrix& A, TVector& x, const TVector& b, const K& w) {
542  TVector v;
543  v=x; //use latest x values in right side calculation
545 
546  }
547  template <typename TVector, typename TMatrix, typename K> //recursion over all matrix rows (A)
548  static void bsorb(const TMatrix& A, TVector& x, TVector& v, const TVector& b, const K& w) {
549  auto rhs = std::get<crow> (b);
550 
551  MultiTypeBlockMatrix_Solver_Col<I,crow,0, TMatrix::M()>::calc_rhs(A,x,v,rhs,w); // calculate right side of equation
552  //solve on blocklevel I-1
553  using M =
554  typename std::remove_cv<
555  typename std::remove_reference<
556  decltype(std::get<crow>( std::get<crow>(A)))
557  >::type
558  >::type;
559  algmeta_itsteps<I-1,M>::bsorb(std::get<crow>( std::get<crow>(A)), std::get<crow>(v),rhs,w);
560  std::get<crow>(x).axpy(w,std::get<crow>(v));
562  }
563 
564 
568  template <typename TVector, typename TMatrix, typename K>
569  static void dbjac(const TMatrix& A, TVector& x, const TVector& b, const K& w) {
570  TVector v(x);
571  v=0; //calc new x in v
573  x.axpy(w,v); //improve x
574  }
575  template <typename TVector, typename TMatrix, typename K>
576  static void dbjac(const TMatrix& A, TVector& x, TVector& v, const TVector& b, const K& w) {
577  auto rhs = std::get<crow> (b);
578 
579  MultiTypeBlockMatrix_Solver_Col<I,crow,0, TMatrix::M()>::calc_rhs(A,x,v,rhs,w); // calculate right side of equation
580  //solve on blocklevel I-1
581  using M =
582  typename std::remove_cv<
583  typename std::remove_reference<
584  decltype(std::get<crow>( std::get<crow>(A)))
585  >::type
586  >::type;
587  algmeta_itsteps<I-1,M>::dbjac(std::get<crow>( std::get<crow>(A)), std::get<crow>(v),rhs,w);
589  }
590 
591 
592 
593 
594  };
595  template<int I, int crow> //recursion end for remain_row = 0
596  class MultiTypeBlockMatrix_Solver<I,crow,0> {
597  public:
598  template <typename TVector, typename TMatrix, typename K>
599  static void dbgs(const TMatrix&, TVector&, TVector&,
600  const TVector&, const K&) {}
601 
602  template <typename TVector, typename TMatrix, typename K>
603  static void bsorf(const TMatrix&, TVector&, TVector&,
604  const TVector&, const K&) {}
605 
606  template <typename TVector, typename TMatrix, typename K>
607  static void bsorb(const TMatrix&, TVector&, TVector&,
608  const TVector&, const K&) {}
609 
610  template <typename TVector, typename TMatrix, typename K>
611  static void dbjac(const TMatrix&, TVector&, TVector&,
612  const TVector&, const K&) {}
613  };
614 
615 } // end namespace Dune
616 
617 namespace std
618 {
623  template <size_t i, typename... Args>
624  struct tuple_element<i,Dune::MultiTypeBlockMatrix<Args...> >
625  {
626  using type = typename std::tuple_element<i, std::tuple<Args...> >::type;
627  };
628 
633  template <typename... Args>
634  struct tuple_size<Dune::MultiTypeBlockMatrix<Args...> >
635  : std::integral_constant<std::size_t, sizeof...(Args)>
636  {};
637 }
638 #endif
part of solvers for MultiTypeBlockVector & MultiTypeBlockMatrix types
Definition: multitypeblockmatrix.hh:450
solver for MultiTypeBlockVector & MultiTypeBlockMatrix types
Definition: multitypeblockmatrix.hh:478
A Matrix class to support different block types.
Definition: multitypeblockmatrix.hh:46
constexpr index_constant< 0 > _0
Compile time index with value 0.
Definition: indices.hh:52
std::integral_constant< std::size_t, i > index_constant
An index constant with value i.
Definition: indices.hh:29
MultiTypeBlockMatrix< FirstRow, Args... > type
Definition: multitypeblockmatrix.hh:56
static void dbjac(const TMatrix &A, TVector &x, const TVector &b, const K &w)
Definition: multitypeblockmatrix.hh:569
void mv(const X &x, Y &y) const
y = A x
Definition: multitypeblockmatrix.hh:185
void mmtv(const X &x, Y &y) const
y -= A^T x
Definition: multitypeblockmatrix.hh:265
void mmhv(const X &x, Y &y) const
y -= A^H x
Definition: multitypeblockmatrix.hh:310
static void dbgs(const TMatrix &A, TVector &x, const TVector &b, const K &w)
Definition: multitypeblockmatrix.hh:485
static void bsorb(const TMatrix &A, TVector &x, const TVector &b, const K &w)
Definition: multitypeblockmatrix.hh:541
MultiTypeBlockMatrix & operator/=(const field_type &k)
vector space division by scalar
Definition: multitypeblockmatrix.hh:141
void usmhv(const field_type &alpha, const X &x, Y &y) const
y += alpha A^H x
Definition: multitypeblockmatrix.hh:325
void usmv(const AlphaType &alpha, const X &x, Y &y) const
y += alpha A x
Definition: multitypeblockmatrix.hh:225
void umhv(const X &x, Y &y) const
y += A^H x
Definition: multitypeblockmatrix.hh:295
static constexpr size_type N()
Return the number of matrix rows.
Definition: multitypeblockmatrix.hh:64
void usmtv(const field_type &alpha, const X &x, Y &y) const
y += alpha A^T x
Definition: multitypeblockmatrix.hh:280
void operator=(const T &newval)
Definition: multitypeblockmatrix.hh:116
void umv(const X &x, Y &y) const
y += A x
Definition: multitypeblockmatrix.hh:195
static void calc_rhs(const TMatrix &A, TVector &x, TVector &v, Trhs &b, const K &w)
Definition: multitypeblockmatrix.hh:456
static void bsorf(const TMatrix &A, TVector &x, const TVector &b, const K &w)
Definition: multitypeblockmatrix.hh:514
auto infinity_norm_real() const
Bastardized version of the infinity-norm / row-sum norm.
Definition: multitypeblockmatrix.hh:385
auto frobenius_norm2() const
square of frobenius norm, need for block recursion
Definition: multitypeblockmatrix.hh:341
std::size_t size_type
Type used for sizes.
Definition: multitypeblockmatrix.hh:59
auto operator[]([[maybe_unused]] const std::integral_constant< size_type, index > indexVariable) -> decltype(std::get< index >(*this))
Random-access operator.
Definition: multitypeblockmatrix.hh:93
auto infinity_norm() const
Bastardized version of the infinity-norm / row-sum norm.
Definition: multitypeblockmatrix.hh:364
FieldTraits< field_type >::real_type frobenius_norm() const
frobenius norm: sqrt(sum over squared values of entries)
Definition: multitypeblockmatrix.hh:358
MultiTypeBlockMatrix & operator+=(const MultiTypeBlockMatrix &b)
Add the entries of another matrix to this one.
Definition: multitypeblockmatrix.hh:157
void mtv(const X &x, Y &y) const
y = A^T x
Definition: multitypeblockmatrix.hh:240
static constexpr size_type M()
Return the number of matrix columns.
Definition: multitypeblockmatrix.hh:70
void mmv(const X &x, Y &y) const
y -= A x
Definition: multitypeblockmatrix.hh:210
void umtv(const X &x, Y &y) const
y += A^T x
Definition: multitypeblockmatrix.hh:250
MultiTypeBlockMatrix & operator*=(const field_type &k)
vector space multiplication with scalar
Definition: multitypeblockmatrix.hh:130
MultiTypeBlockMatrix & operator-=(const MultiTypeBlockMatrix &b)
Subtract the entries of another matrix from this one.
Definition: multitypeblockmatrix.hh:172
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
void bsorb(const M &A, X &x, const Y &b, const K &w)
SSOR step.
Definition: gsetc.hh:646
void dbjac(const M &A, X &x, const Y &b, const K &w)
Jacobi step.
Definition: gsetc.hh:658
void dbgs(const M &A, X &x, const Y &b, const K &w)
GS step.
Definition: gsetc.hh:622
void bsorf(const M &A, X &x, const Y &b, const K &w)
SOR step.
Definition: gsetc.hh:634
Simple iterative methods like Jacobi, Gauss-Seidel, SOR, SSOR, etc. in a generic way.
Dune namespace.
Definition: alignedallocator.hh:13
constexpr std::integral_constant< std::size_t, sizeof...(II)> size(std::integer_sequence< T, II... >)
Return the size of the sequence.
Definition: integersequence.hh:75
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.80.0 (Apr 27, 22:29, 2024)