Dune Core Modules (2.9.0)

basearray.hh
Go to the documentation of this file.
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_BASEARRAY_HH
6 #define DUNE_ISTL_BASEARRAY_HH
7 
8 #include "assert.h"
9 #include <cmath>
10 #include <cstddef>
11 #include <memory>
12 #include <algorithm>
13 
14 #include "istlexception.hh"
16 
21 namespace Dune {
22 
24 namespace Imp {
25 
50  template<class B, class A=std::allocator<B> >
51  class base_array_unmanaged
52  {
53  public:
54 
55  //===== type definitions and constants
56 
58  typedef B member_type;
59 
61  typedef A allocator_type;
62 
64  typedef typename A::size_type size_type;
65 
67  using reference = B&;
68 
70  using const_reference = const B&;
71 
72  //===== access to components
73 
75  reference operator[] (size_type i)
76  {
77 #ifdef DUNE_ISTL_WITH_CHECKING
78  if (i>=n) DUNE_THROW(ISTLError,"index out of range");
79 #endif
80  return p[i];
81  }
82 
84  const_reference operator[] (size_type i) const
85  {
86 #ifdef DUNE_ISTL_WITH_CHECKING
87  if (i>=n) DUNE_THROW(ISTLError,"index out of range");
88 #endif
89  return p[i];
90  }
91 
93  template<class T>
94  class RealIterator
95  : public RandomAccessIteratorFacade<RealIterator<T>, T>
96  {
97  public:
99  typedef typename std::remove_const<T>::type ValueType;
100 
101  friend class RandomAccessIteratorFacade<RealIterator<const ValueType>, const ValueType>;
102  friend class RandomAccessIteratorFacade<RealIterator<ValueType>, ValueType>;
103  friend class RealIterator<const ValueType>;
104  friend class RealIterator<ValueType>;
105 
107  RealIterator ()
108  : p(0), i(0)
109  {}
110 
111  RealIterator (const B* _p, B* _i) : p(_p), i(_i)
112  { }
113 
114  RealIterator(const RealIterator<ValueType>& it)
115  : p(it.p), i(it.i)
116  {}
117 
119  size_type index () const
120  {
121  return i-p;
122  }
123 
125  bool equals (const RealIterator<ValueType>& other) const
126  {
127  assert(other.p==p);
128  return i==other.i;
129  }
130 
132  bool equals (const RealIterator<const ValueType>& other) const
133  {
134  assert(other.p==p);
135  return i==other.i;
136  }
137 
138  std::ptrdiff_t distanceTo(const RealIterator& o) const
139  {
140  return o.i-i;
141  }
142 
143  private:
145  void increment()
146  {
147  ++i;
148  }
149 
151  void decrement()
152  {
153  --i;
154  }
155 
156  // Needed for operator[] of the iterator
157  reference elementAt (std::ptrdiff_t offset) const
158  {
159  return *(i+offset);
160  }
161 
163  reference dereference () const
164  {
165  return *i;
166  }
167 
168  void advance(std::ptrdiff_t d)
169  {
170  i+=d;
171  }
172 
173  const B* p;
174  B* i;
175  };
176 
178  typedef RealIterator<B> iterator;
179 
180 
182  iterator begin ()
183  {
184  return iterator(p,p);
185  }
186 
188  iterator end ()
189  {
190  return iterator(p,p+n);
191  }
192 
195  iterator beforeEnd ()
196  {
197  return iterator(p,p+n-1);
198  }
199 
202  iterator beforeBegin ()
203  {
204  return iterator(p,p-1);
205  }
206 
208  iterator find (size_type i)
209  {
210  return iterator(p,p+std::min(i,n));
211  }
212 
214  typedef RealIterator<const B> const_iterator;
215 
217  const_iterator begin () const
218  {
219  return const_iterator(p,p+0);
220  }
221 
223  const_iterator end () const
224  {
225  return const_iterator(p,p+n);
226  }
227 
230  const_iterator beforeEnd () const
231  {
232  return const_iterator(p,p+n-1);
233  }
234 
237  const_iterator beforeBegin () const
238  {
239  return const_iterator(p,p-1);
240  }
241 
243  const_iterator find (size_type i) const
244  {
245  return const_iterator(p,p+std::min(i,n));
246  }
247 
248 
249  //===== sizes
250 
252  size_type size () const
253  {
254  return n;
255  }
256 
258  const B* data() const
259  {
260  return p;
261  }
262 
264  B* data()
265  {
266  return p;
267  }
268 
269  protected:
271  base_array_unmanaged ()
272  : n(0), p(0)
273  {}
275  base_array_unmanaged (size_type n_, B* p_)
276  : n(n_), p(p_)
277  {}
278  size_type n; // number of elements in array
279  B *p; // pointer to dynamically allocated built-in array
280  };
281 
282 
283 
305  template<class B, class A=std::allocator<B> >
306  class compressed_base_array_unmanaged
307  {
308  public:
309 
310  //===== type definitions and constants
311 
313  typedef B member_type;
314 
316  typedef A allocator_type;
317 
319  typedef typename A::size_type size_type;
320 
322  using reference = B&;
323 
325  using const_reference = const B&;
326 
327  //===== access to components
328 
330  reference operator[] (size_type i)
331  {
332  const size_type* lb = std::lower_bound(j, j+n, i);
333  if (lb == j+n || *lb != i)
334  DUNE_THROW(ISTLError,"index "<<i<<" not in compressed array");
335  return p[lb-j];
336  }
337 
339  const_reference operator[] (size_type i) const
340  {
341  const size_type* lb = std::lower_bound(j, j+n, i);
342  if (lb == j+n || *lb != i)
343  DUNE_THROW(ISTLError,"index "<<i<<" not in compressed array");
344  return p[lb-j];
345  }
346 
348  template<class T>
349  class RealIterator
350  : public BidirectionalIteratorFacade<RealIterator<T>, T>
351  {
352  public:
354  typedef typename std::remove_const<T>::type ValueType;
355 
356  friend class BidirectionalIteratorFacade<RealIterator<const ValueType>, const ValueType>;
357  friend class BidirectionalIteratorFacade<RealIterator<ValueType>, ValueType>;
358  friend class RealIterator<const ValueType>;
359  friend class RealIterator<ValueType>;
360 
362  RealIterator ()
363  : p(0), j(0), i(0)
364  {}
365 
367  RealIterator (B* _p, size_type* _j, size_type _i)
368  : p(_p), j(_j), i(_i)
369  { }
370 
374  RealIterator(const RealIterator<ValueType>& it)
375  : p(it.p), j(it.j), i(it.i)
376  {}
377 
378 
380  bool equals (const RealIterator<ValueType>& it) const
381  {
382  assert(p==it.p);
383  return (i)==(it.i);
384  }
385 
387  bool equals (const RealIterator<const ValueType>& it) const
388  {
389  assert(p==it.p);
390  return (i)==(it.i);
391  }
392 
393 
395  size_type index () const
396  {
397  return j[i];
398  }
399 
401  void setindex (size_type k)
402  {
403  return j[i] = k;
404  }
405 
413  size_type offset () const
414  {
415  return i;
416  }
417 
418  private:
420  void increment()
421  {
422  ++i;
423  }
424 
426  void decrement()
427  {
428  --i;
429  }
430 
432  reference dereference () const
433  {
434  return p[i];
435  }
436 
437  B* p;
438  size_type* j;
439  size_type i;
440  };
441 
443  typedef RealIterator<B> iterator;
444 
446  iterator begin ()
447  {
448  return iterator(p,j,0);
449  }
450 
452  iterator end ()
453  {
454  return iterator(p,j,n);
455  }
456 
459  iterator beforeEnd ()
460  {
461  return iterator(p,j,n-1);
462  }
463 
466  iterator beforeBegin ()
467  {
468  return iterator(p,j,-1);
469  }
470 
472  iterator find (size_type i)
473  {
474  const size_type* lb = std::lower_bound(j, j+n, i);
475  return (lb != j+n && *lb == i)
476  ? iterator(p,j,lb-j)
477  : end();
478  }
479 
481  typedef RealIterator<const B> const_iterator;
482 
484  const_iterator begin () const
485  {
486  return const_iterator(p,j,0);
487  }
488 
490  const_iterator end () const
491  {
492  return const_iterator(p,j,n);
493  }
494 
497  const_iterator beforeEnd () const
498  {
499  return const_iterator(p,j,n-1);
500  }
501 
504  const_iterator beforeBegin () const
505  {
506  return const_iterator(p,j,-1);
507  }
508 
510  const_iterator find (size_type i) const
511  {
512  const size_type* lb = std::lower_bound(j, j+n, i);
513  return (lb != j+n && *lb == i)
514  ? const_iterator(p,j,lb-j)
515  : end();
516  }
517 
518  //===== sizes
519 
521  size_type size () const
522  {
523  return n;
524  }
525 
526  protected:
528  compressed_base_array_unmanaged ()
529  : n(0), p(0), j(0)
530  {}
531 
532  size_type n; // number of elements in array
533  B *p; // pointer to dynamically allocated built-in array
534  size_type* j; // the index set
535  };
536 
537 } // end namespace Imp
538 
539 } // end namespace
540 
541 #endif
#define DUNE_THROW(E, m)
Definition: exceptions.hh:218
constexpr auto equals(T1 &&t1, T2 &&t2)
Equality comparison.
Definition: hybridutilities.hh:402
constexpr decltype(auto) elementAt(Container &&c, Index &&i)
Get element at given position from container.
Definition: hybridutilities.hh:135
auto min(ADLTag< 0 >, const V &v1, const V &v2)
implements binary Simd::min()
Definition: defaults.hh:89
This file implements iterator facade classes for writing stl conformant iterators.
Dune namespace.
Definition: alignedallocator.hh:13
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.80.0 (May 14, 22:30, 2024)