Dune Core Modules (2.7.1)

basearray.hh
Go to the documentation of this file.
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_ISTL_BASEARRAY_HH
4 #define DUNE_ISTL_BASEARRAY_HH
5 
6 #include "assert.h"
7 #include <cmath>
8 #include <cstddef>
9 #include <memory>
10 #include <algorithm>
11 
12 #include "istlexception.hh"
14 
19 namespace Dune {
20 
22 namespace Imp {
23 
48  template<class B, class A=std::allocator<B> >
49  class base_array_unmanaged
50  {
51  public:
52 
53  //===== type definitions and constants
54 
56  typedef B member_type;
57 
59  typedef A allocator_type;
60 
62  typedef typename A::size_type size_type;
63 
65  using reference = B&;
66 
68  using const_reference = const B&;
69 
70  //===== access to components
71 
73  reference operator[] (size_type i)
74  {
75 #ifdef DUNE_ISTL_WITH_CHECKING
76  if (i>=n) DUNE_THROW(ISTLError,"index out of range");
77 #endif
78  return p[i];
79  }
80 
82  const_reference operator[] (size_type i) const
83  {
84 #ifdef DUNE_ISTL_WITH_CHECKING
85  if (i>=n) DUNE_THROW(ISTLError,"index out of range");
86 #endif
87  return p[i];
88  }
89 
91  template<class T>
92  class RealIterator
93  : public RandomAccessIteratorFacade<RealIterator<T>, T>
94  {
95  public:
97  typedef typename std::remove_const<T>::type ValueType;
98 
99  friend class RandomAccessIteratorFacade<RealIterator<const ValueType>, const ValueType>;
100  friend class RandomAccessIteratorFacade<RealIterator<ValueType>, ValueType>;
101  friend class RealIterator<const ValueType>;
102  friend class RealIterator<ValueType>;
103 
105  RealIterator ()
106  : p(0), i(0)
107  {}
108 
109  RealIterator (const B* _p, B* _i) : p(_p), i(_i)
110  { }
111 
112  RealIterator(const RealIterator<ValueType>& it)
113  : p(it.p), i(it.i)
114  {}
115 
117  size_type index () const
118  {
119  return i-p;
120  }
121 
123  bool equals (const RealIterator<ValueType>& other) const
124  {
125  assert(other.p==p);
126  return i==other.i;
127  }
128 
130  bool equals (const RealIterator<const ValueType>& other) const
131  {
132  assert(other.p==p);
133  return i==other.i;
134  }
135 
136  std::ptrdiff_t distanceTo(const RealIterator& o) const
137  {
138  return o.i-i;
139  }
140 
141  private:
143  void increment()
144  {
145  ++i;
146  }
147 
149  void decrement()
150  {
151  --i;
152  }
153 
154  // Needed for operator[] of the iterator
155  reference elementAt (std::ptrdiff_t offset) const
156  {
157  return *(i+offset);
158  }
159 
161  reference dereference () const
162  {
163  return *i;
164  }
165 
166  void advance(std::ptrdiff_t d)
167  {
168  i+=d;
169  }
170 
171  const B* p;
172  B* i;
173  };
174 
176  typedef RealIterator<B> iterator;
177 
178 
180  iterator begin ()
181  {
182  return iterator(p,p);
183  }
184 
186  iterator end ()
187  {
188  return iterator(p,p+n);
189  }
190 
193  iterator beforeEnd ()
194  {
195  return iterator(p,p+n-1);
196  }
197 
200  iterator beforeBegin ()
201  {
202  return iterator(p,p-1);
203  }
204 
206  iterator find (size_type i)
207  {
208  return iterator(p,p+std::min(i,n));
209  }
210 
212  typedef RealIterator<const B> const_iterator;
213 
215  const_iterator begin () const
216  {
217  return const_iterator(p,p+0);
218  }
219 
221  const_iterator end () const
222  {
223  return const_iterator(p,p+n);
224  }
225 
228  const_iterator beforeEnd () const
229  {
230  return const_iterator(p,p+n-1);
231  }
232 
235  const_iterator beforeBegin () const
236  {
237  return const_iterator(p,p-1);
238  }
239 
241  const_iterator find (size_type i) const
242  {
243  return const_iterator(p,p+std::min(i,n));
244  }
245 
246 
247  //===== sizes
248 
250  size_type size () const
251  {
252  return n;
253  }
254 
256  const B* data() const
257  {
258  return p;
259  }
260 
262  B* data()
263  {
264  return p;
265  }
266 
267  protected:
269  base_array_unmanaged ()
270  : n(0), p(0)
271  {}
273  base_array_unmanaged (size_type n_, B* p_)
274  : n(n_), p(p_)
275  {}
276  size_type n; // number of elements in array
277  B *p; // pointer to dynamically allocated built-in array
278  };
279 
280 
281 
303  template<class B, class A=std::allocator<B> >
304  class compressed_base_array_unmanaged
305  {
306  public:
307 
308  //===== type definitions and constants
309 
311  typedef B member_type;
312 
314  typedef A allocator_type;
315 
317  typedef typename A::size_type size_type;
318 
320  using reference = B&;
321 
323  using const_reference = const B&;
324 
325  //===== access to components
326 
328  reference operator[] (size_type i)
329  {
330  const size_type* lb = std::lower_bound(j, j+n, i);
331  if (lb == j+n || *lb != i)
332  DUNE_THROW(ISTLError,"index "<<i<<" not in compressed array");
333  return p[lb-j];
334  }
335 
337  const_reference operator[] (size_type i) const
338  {
339  const size_type* lb = std::lower_bound(j, j+n, i);
340  if (lb == j+n || *lb != i)
341  DUNE_THROW(ISTLError,"index "<<i<<" not in compressed array");
342  return p[lb-j];
343  }
344 
346  template<class T>
347  class RealIterator
348  : public BidirectionalIteratorFacade<RealIterator<T>, T>
349  {
350  public:
352  typedef typename std::remove_const<T>::type ValueType;
353 
354  friend class BidirectionalIteratorFacade<RealIterator<const ValueType>, const ValueType>;
355  friend class BidirectionalIteratorFacade<RealIterator<ValueType>, ValueType>;
356  friend class RealIterator<const ValueType>;
357  friend class RealIterator<ValueType>;
358 
360  RealIterator ()
361  : p(0), j(0), i(0)
362  {}
363 
365  RealIterator (B* _p, size_type* _j, size_type _i)
366  : p(_p), j(_j), i(_i)
367  { }
368 
372  RealIterator(const RealIterator<ValueType>& it)
373  : p(it.p), j(it.j), i(it.i)
374  {}
375 
376 
378  bool equals (const RealIterator<ValueType>& it) const
379  {
380  assert(p==it.p);
381  return (i)==(it.i);
382  }
383 
385  bool equals (const RealIterator<const ValueType>& it) const
386  {
387  assert(p==it.p);
388  return (i)==(it.i);
389  }
390 
391 
393  size_type index () const
394  {
395  return j[i];
396  }
397 
399  void setindex (size_type k)
400  {
401  return j[i] = k;
402  }
403 
411  size_type offset () const
412  {
413  return i;
414  }
415 
416  private:
418  void increment()
419  {
420  ++i;
421  }
422 
424  void decrement()
425  {
426  --i;
427  }
428 
430  reference dereference () const
431  {
432  return p[i];
433  }
434 
435  B* p;
436  size_type* j;
437  size_type i;
438  };
439 
441  typedef RealIterator<B> iterator;
442 
444  iterator begin ()
445  {
446  return iterator(p,j,0);
447  }
448 
450  iterator end ()
451  {
452  return iterator(p,j,n);
453  }
454 
457  iterator beforeEnd ()
458  {
459  return iterator(p,j,n-1);
460  }
461 
464  iterator beforeBegin ()
465  {
466  return iterator(p,j,-1);
467  }
468 
470  iterator find (size_type i)
471  {
472  const size_type* lb = std::lower_bound(j, j+n, i);
473  return (lb != j+n && *lb == i)
474  ? iterator(p,j,lb-j)
475  : end();
476  }
477 
479  typedef RealIterator<const B> const_iterator;
480 
482  const_iterator begin () const
483  {
484  return const_iterator(p,j,0);
485  }
486 
488  const_iterator end () const
489  {
490  return const_iterator(p,j,n);
491  }
492 
495  const_iterator beforeEnd () const
496  {
497  return const_iterator(p,j,n-1);
498  }
499 
502  const_iterator beforeBegin () const
503  {
504  return const_iterator(p,j,-1);
505  }
506 
508  const_iterator find (size_type i) const
509  {
510  const size_type* lb = std::lower_bound(j, j+n, i);
511  return (lb != j+n && *lb == i)
512  ? const_iterator(p,j,lb-j)
513  : end();
514  }
515 
516  //===== sizes
517 
519  size_type size () const
520  {
521  return n;
522  }
523 
524  protected:
526  compressed_base_array_unmanaged ()
527  : n(0), p(0), j(0)
528  {}
529 
530  size_type n; // number of elements in array
531  B *p; // pointer to dynamically allocated built-in array
532  size_type* j; // the index set
533  };
534 
535 } // end namespace Imp
536 
537 } // end namespace
538 
539 #endif
#define DUNE_THROW(E, m)
Definition: exceptions.hh:216
constexpr auto equals(T1 &&t1, T2 &&t2)
Equality comparison.
Definition: hybridutilities.hh:401
constexpr decltype(auto) elementAt(Container &&c, Index &&i)
Get element at given position from container.
Definition: hybridutilities.hh:134
auto min(ADLTag< 0 >, const V &v1, const V &v2)
implements binary Simd::min()
Definition: defaults.hh:87
This file implements iterator facade classes for writing stl conformant iterators.
Dune namespace.
Definition: alignedallocator.hh:14
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.80.0 (May 16, 22:29, 2024)