Dune Core Modules (2.6.0)

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 
297  template<class B, class A=std::allocator<B> >
298  class base_array_window : public base_array_unmanaged<B,A>
299  {
300  public:
301 
302  //===== type definitions and constants
303 
305  typedef B member_type;
306 
308  typedef A allocator_type;
309 
311  typedef typename base_array_unmanaged<B,A>::iterator iterator;
312 
314  typedef typename base_array_unmanaged<B,A>::const_iterator const_iterator;
315 
317  typedef typename base_array_unmanaged<B,A>::size_type size_type;
318 
320  typedef typename A::difference_type difference_type;
321 
322  //===== constructors and such
323 
325  base_array_window ()
326  : base_array_unmanaged<B,A>()
327  { }
328 
330  base_array_window (B* _p, size_type _n)
331  : base_array_unmanaged<B,A>(_n ,_p)
332  {}
333 
334  //===== window manipulation methods
335 
337  void set (size_type _n, B* _p)
338  {
339  this->n = _n;
340  this->p = _p;
341  }
342 
344  void advance (difference_type newsize)
345  {
346  this->p += this->n;
347  this->n = newsize;
348  }
349 
351  void move (difference_type offset, size_type newsize)
352  {
353  this->p += offset;
354  this->n = newsize;
355  }
356 
358  void move (difference_type offset)
359  {
360  this->p += offset;
361  }
362 
364  B* getptr ()
365  {
366  return this->p;
367  }
368  };
369 
370 
371 
388  template<class B, class A=std::allocator<B> >
389  class base_array : public base_array_unmanaged<B,A>
390  {
391  public:
392 
393  //===== type definitions and constants
394 
396  typedef B member_type;
397 
399  typedef A allocator_type;
400 
402  typedef typename base_array_unmanaged<B,A>::iterator iterator;
403 
405  typedef typename base_array_unmanaged<B,A>::const_iterator const_iterator;
406 
408  typedef typename base_array_unmanaged<B,A>::size_type size_type;
409 
411  typedef typename A::difference_type difference_type;
412 
413  //===== constructors and such
414 
416  base_array ()
417  : base_array_unmanaged<B,A>()
418  {}
419 
421  base_array (size_type _n)
422  : base_array_unmanaged<B,A>(_n, 0)
423  {
424  if (this->n>0) {
425  this->p = allocator_.allocate(this->n);
426  new (this->p)B[this->n];
427  } else
428  {
429  this->n = 0;
430  this->p = 0;
431  }
432  }
433 
435  base_array (const base_array& a)
436  {
437  // allocate memory with same size as a
438  this->n = a.n;
439 
440  if (this->n>0) {
441  this->p = allocator_.allocate(this->n);
442  new (this->p)B[this->n];
443  } else
444  {
445  this->n = 0;
446  this->p = 0;
447  }
448 
449  // and copy elements
450  for (size_type i=0; i<this->n; i++) this->p[i]=a.p[i];
451  }
452 
454  ~base_array ()
455  {
456  if (this->n>0) {
457  int i=this->n;
458  while (i)
459  this->p[--i].~B();
460  allocator_.deallocate(this->p,this->n);
461  }
462  }
463 
465  void resize (size_type _n)
466  {
467  if (this->n==_n) return;
468 
469  if (this->n>0) {
470  int i=this->n;
471  while (i)
472  this->p[--i].~B();
473  allocator_.deallocate(this->p,this->n);
474  }
475  this->n = _n;
476  if (this->n>0) {
477  this->p = allocator_.allocate(this->n);
478  new (this->p)B[this->n];
479  } else
480  {
481  this->n = 0;
482  this->p = 0;
483  }
484  }
485 
487  base_array& operator= (const base_array& a)
488  {
489  if (&a!=this) // check if this and a are different objects
490  {
491  // adjust size of array
492  if (this->n!=a.n) // check if size is different
493  {
494  if (this->n>0) {
495  int i=this->n;
496  while (i)
497  this->p[--i].~B();
498  allocator_.deallocate(this->p,this->n); // delete old memory
499  }
500  this->n = a.n;
501  if (this->n>0) {
502  this->p = allocator_.allocate(this->n);
503  new (this->p)B[this->n];
504  } else
505  {
506  this->n = 0;
507  this->p = 0;
508  }
509  }
510  // copy data
511  for (size_type i=0; i<this->n; i++) this->p[i]=a.p[i];
512  }
513  return *this;
514  }
515 
516  protected:
517 
518  A allocator_;
519  };
520 
521 
522 
523 
545  template<class B, class A=std::allocator<B> >
546  class compressed_base_array_unmanaged
547  {
548  public:
549 
550  //===== type definitions and constants
551 
553  typedef B member_type;
554 
556  typedef A allocator_type;
557 
559  typedef typename A::size_type size_type;
560 
562  using reference = B&;
563 
565  using const_reference = const B&;
566 
567  //===== access to components
568 
570  reference operator[] (size_type i)
571  {
572  const size_type* lb = std::lower_bound(j, j+n, i);
573  if (lb == j+n || *lb != i)
574  DUNE_THROW(ISTLError,"index "<<i<<" not in compressed array");
575  return p[lb-j];
576  }
577 
579  const_reference operator[] (size_type i) const
580  {
581  const size_type* lb = std::lower_bound(j, j+n, i);
582  if (lb == j+n || *lb != i)
583  DUNE_THROW(ISTLError,"index "<<i<<" not in compressed array");
584  return p[lb-j];
585  }
586 
588  template<class T>
589  class RealIterator
590  : public BidirectionalIteratorFacade<RealIterator<T>, T>
591  {
592  public:
594  typedef typename std::remove_const<T>::type ValueType;
595 
596  friend class BidirectionalIteratorFacade<RealIterator<const ValueType>, const ValueType>;
597  friend class BidirectionalIteratorFacade<RealIterator<ValueType>, ValueType>;
598  friend class RealIterator<const ValueType>;
599  friend class RealIterator<ValueType>;
600 
602  RealIterator ()
603  : p(0), j(0), i(0)
604  {}
605 
607  RealIterator (B* _p, size_type* _j, size_type _i)
608  : p(_p), j(_j), i(_i)
609  { }
610 
614  RealIterator(const RealIterator<ValueType>& it)
615  : p(it.p), j(it.j), i(it.i)
616  {}
617 
618 
620  bool equals (const RealIterator<ValueType>& it) const
621  {
622  assert(p==it.p);
623  return (i)==(it.i);
624  }
625 
627  bool equals (const RealIterator<const ValueType>& it) const
628  {
629  assert(p==it.p);
630  return (i)==(it.i);
631  }
632 
633 
635  size_type index () const
636  {
637  return j[i];
638  }
639 
641  void setindex (size_type k)
642  {
643  return j[i] = k;
644  }
645 
653  size_type offset () const
654  {
655  return i;
656  }
657 
658  private:
660  void increment()
661  {
662  ++i;
663  }
664 
666  void decrement()
667  {
668  --i;
669  }
670 
672  reference dereference () const
673  {
674  return p[i];
675  }
676 
677  B* p;
678  size_type* j;
679  size_type i;
680  };
681 
683  typedef RealIterator<B> iterator;
684 
686  iterator begin ()
687  {
688  return iterator(p,j,0);
689  }
690 
692  iterator end ()
693  {
694  return iterator(p,j,n);
695  }
696 
699  iterator beforeEnd ()
700  {
701  return iterator(p,j,n-1);
702  }
703 
706  iterator beforeBegin ()
707  {
708  return iterator(p,j,-1);
709  }
710 
712  iterator find (size_type i)
713  {
714  const size_type* lb = std::lower_bound(j, j+n, i);
715  return (lb != j+n && *lb == i)
716  ? iterator(p,j,lb-j)
717  : end();
718  }
719 
721  typedef RealIterator<const B> const_iterator;
722 
724  const_iterator begin () const
725  {
726  return const_iterator(p,j,0);
727  }
728 
730  const_iterator end () const
731  {
732  return const_iterator(p,j,n);
733  }
734 
737  const_iterator beforeEnd () const
738  {
739  return const_iterator(p,j,n-1);
740  }
741 
744  const_iterator beforeBegin () const
745  {
746  return const_iterator(p,j,-1);
747  }
748 
750  const_iterator find (size_type i) const
751  {
752  const size_type* lb = std::lower_bound(j, j+n, i);
753  return (lb != j+n && *lb == i)
754  ? const_iterator(p,j,lb-j)
755  : end();
756  }
757 
758  //===== sizes
759 
761  size_type size () const
762  {
763  return n;
764  }
765 
766  protected:
768  compressed_base_array_unmanaged ()
769  : n(0), p(0), j(0)
770  {}
771 
772  size_type n; // number of elements in array
773  B *p; // pointer to dynamically allocated built-in array
774  size_type* j; // the index set
775  };
776 
777 } // end namespace Imp
778 
779 } // end namespace
780 
781 #endif
#define DUNE_THROW(E, m)
Definition: exceptions.hh:216
constexpr auto equals(T1 &&t1, T2 &&t2)
Equality comparison.
Definition: hybridutilities.hh:435
constexpr decltype(auto) elementAt(Container &&c, Index &&i)
Get element at given position from container.
Definition: hybridutilities.hh:133
This file implements iterator facade classes for writing stl conformant iterators.
Dune namespace.
Definition: alignedallocator.hh:10
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.80.0 (May 5, 22:29, 2024)