Dune Core Modules (unstable)

basearray.hh
Go to the documentation of this file.
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_BASEARRAY_HH
6#define DUNE_ISTL_BASEARRAY_HH
7
8#include <cassert>
9#include <cmath>
10#include <cstddef>
11#include <memory>
12#include <algorithm>
13#include <type_traits>
14
15#include "istlexception.hh"
17
22namespace Dune {
23
25namespace Imp {
26
47 template<class B, class ST=std::size_t >
48 class base_array_unmanaged
49 {
50 public:
51
52 //===== type definitions and constants
53
55 typedef B member_type;
56
58 typedef ST size_type;
59
61 using reference = B&;
62
64 using const_reference = const B&;
65
66 //===== access to components
67
69 reference operator[] (size_type i)
70 {
71#ifdef DUNE_ISTL_WITH_CHECKING
72 if (i>=n) DUNE_THROW(ISTLError,"index out of range");
73#endif
74 return p[i];
75 }
76
78 const_reference operator[] (size_type i) const
79 {
80#ifdef DUNE_ISTL_WITH_CHECKING
81 if (i>=n) DUNE_THROW(ISTLError,"index out of range");
82#endif
83 return p[i];
84 }
85
87 template<class T>
88 class RealIterator
89 : public RandomAccessIteratorFacade<RealIterator<T>, T>
90 {
91 public:
93 typedef typename std::remove_const<T>::type ValueType;
94
95 friend class RandomAccessIteratorFacade<RealIterator<const ValueType>, const ValueType>;
96 friend class RandomAccessIteratorFacade<RealIterator<ValueType>, ValueType>;
97 friend class RealIterator<const ValueType>;
98 friend class RealIterator<ValueType>;
99
101 RealIterator () = default;
102
103 RealIterator (const B* _p, B* _i)
104 : p(_p), i(_i)
105 {}
106
107 template <class T_,
108 std::enable_if_t<std::is_same_v<std::remove_const_t<T>, std::remove_const_t<T_>>, int> = 0>
109 RealIterator (const RealIterator<T_>& other)
110 : p(other.p), i(other.i)
111 {}
112
113 template <class T_,
114 std::enable_if_t<std::is_same_v<std::remove_const_t<T>, std::remove_const_t<T_>>, int> = 0>
115 RealIterator& operator= (const RealIterator<T_>& other)
116 {
117 p = other.p;
118 i = other.i;
119 return *this;
120 }
121
123 size_type index () const
124 {
125 return i-p;
126 }
127
129 bool equals (const RealIterator<ValueType>& other) const
130 {
131 assert(other.p==p);
132 return i==other.i;
133 }
134
136 bool equals (const RealIterator<const ValueType>& other) const
137 {
138 assert(other.p==p);
139 return i==other.i;
140 }
141
142 std::ptrdiff_t distanceTo(const RealIterator& o) const
143 {
144 return o.i-i;
145 }
146
147 private:
149 void increment()
150 {
151 ++i;
152 }
153
155 void decrement()
156 {
157 --i;
158 }
159
160 // Needed for operator[] of the iterator
161 reference elementAt (std::ptrdiff_t offset) const
162 {
163 return *(i+offset);
164 }
165
167 reference dereference () const
168 {
169 return *i;
170 }
171
172 void advance(std::ptrdiff_t d)
173 {
174 i+=d;
175 }
176
177 const B* p = nullptr;
178 B* i = nullptr;
179 };
180
182 typedef RealIterator<B> iterator;
183
184
186 iterator begin ()
187 {
188 return iterator(p,p);
189 }
190
192 iterator end ()
193 {
194 return iterator(p,p+n);
195 }
196
199 iterator beforeEnd ()
200 {
201 return iterator(p,p+n-1);
202 }
203
206 iterator beforeBegin ()
207 {
208 return iterator(p,p-1);
209 }
210
212 iterator find (size_type i)
213 {
214 return iterator(p,p+std::min(i,n));
215 }
216
218 typedef RealIterator<const B> const_iterator;
219
221 const_iterator begin () const
222 {
223 return const_iterator(p,p+0);
224 }
225
227 const_iterator end () const
228 {
229 return const_iterator(p,p+n);
230 }
231
234 const_iterator beforeEnd () const
235 {
236 return const_iterator(p,p+n-1);
237 }
238
241 const_iterator beforeBegin () const
242 {
243 return const_iterator(p,p-1);
244 }
245
247 const_iterator find (size_type i) const
248 {
249 return const_iterator(p,p+std::min(i,n));
250 }
251
252
253 //===== sizes
254
256 size_type size () const
257 {
258 return n;
259 }
260
262 const B* data() const
263 {
264 return p;
265 }
266
268 B* data()
269 {
270 return p;
271 }
272
273 protected:
275 base_array_unmanaged ()
276 : n(0), p(0)
277 {}
279 base_array_unmanaged (size_type n_, B* p_)
280 : n(n_), p(p_)
281 {}
282 size_type n; // number of elements in array
283 B *p; // pointer to dynamically allocated built-in array
284 };
285
286
287
309 template<class B, class ST=std::size_t >
310 class compressed_base_array_unmanaged
311 {
312 public:
313
314 //===== type definitions and constants
315
317 typedef B member_type;
318
320 typedef ST size_type;
321
323 using reference = B&;
324
326 using const_reference = const B&;
327
328 //===== access to components
329
331 reference operator[] (size_type i)
332 {
333 const size_type* lb = std::lower_bound(j, j+n, i);
334 if (lb == j+n || *lb != i)
335 DUNE_THROW(ISTLError,"index "<<i<<" not in compressed array");
336 return p[lb-j];
337 }
338
340 const_reference operator[] (size_type i) const
341 {
342 const size_type* lb = std::lower_bound(j, j+n, i);
343 if (lb == j+n || *lb != i)
344 DUNE_THROW(ISTLError,"index "<<i<<" not in compressed array");
345 return p[lb-j];
346 }
347
349 template<class T>
350 class RealIterator
351 : public BidirectionalIteratorFacade<RealIterator<T>, T>
352 {
353 public:
355 typedef typename std::remove_const<T>::type ValueType;
356
357 friend class BidirectionalIteratorFacade<RealIterator<const ValueType>, const ValueType>;
358 friend class BidirectionalIteratorFacade<RealIterator<ValueType>, ValueType>;
359 friend class RealIterator<const ValueType>;
360 friend class RealIterator<ValueType>;
361
363 RealIterator () = default;
364
366 RealIterator (B* _p, size_type* _j, size_type _i)
367 : p(_p), j(_j), i(_i)
368 {}
369
370 template <class T_,
371 std::enable_if_t<std::is_same_v<std::remove_const_t<T>, std::remove_const_t<T_>>, int> = 0>
372 RealIterator (const RealIterator<T_>& other)
373 : p(other.p), j(other.j), i(other.i)
374 {}
375
376 template <class T_,
377 std::enable_if_t<std::is_same_v<std::remove_const_t<T>, std::remove_const_t<T_>>, int> = 0>
378 RealIterator& operator= (const RealIterator<T_>& other)
379 {
380 p = other.p;
381 j = other.j;
382 i = other.i;
383 return *this;
384 }
385
386
388 bool equals (const RealIterator<ValueType>& it) const
389 {
390 assert(p==it.p);
391 return (i)==(it.i);
392 }
393
395 bool equals (const RealIterator<const ValueType>& it) const
396 {
397 assert(p==it.p);
398 return (i)==(it.i);
399 }
400
401
403 size_type index () const
404 {
405 return j[i];
406 }
407
409 void setindex (size_type k)
410 {
411 return j[i] = k;
412 }
413
421 size_type offset () const
422 {
423 return i;
424 }
425
426 private:
428 void increment()
429 {
430 ++i;
431 }
432
434 void decrement()
435 {
436 --i;
437 }
438
440 reference dereference () const
441 {
442 return p[i];
443 }
444
445 B* p = nullptr;
446 size_type* j = nullptr;
447 size_type i = 0;
448 };
449
451 typedef RealIterator<B> iterator;
452
454 iterator begin ()
455 {
456 return iterator(p,j,0);
457 }
458
460 iterator end ()
461 {
462 return iterator(p,j,n);
463 }
464
467 iterator beforeEnd ()
468 {
469 return iterator(p,j,n-1);
470 }
471
474 iterator beforeBegin ()
475 {
476 return iterator(p,j,-1);
477 }
478
480 iterator find (size_type i)
481 {
482 const size_type* lb = std::lower_bound(j, j+n, i);
483 return (lb != j+n && *lb == i)
484 ? iterator(p,j,lb-j)
485 : end();
486 }
487
489 typedef RealIterator<const B> const_iterator;
490
492 const_iterator begin () const
493 {
494 return const_iterator(p,j,0);
495 }
496
498 const_iterator end () const
499 {
500 return const_iterator(p,j,n);
501 }
502
505 const_iterator beforeEnd () const
506 {
507 return const_iterator(p,j,n-1);
508 }
509
512 const_iterator beforeBegin () const
513 {
514 return const_iterator(p,j,-1);
515 }
516
518 const_iterator find (size_type i) const
519 {
520 const size_type* lb = std::lower_bound(j, j+n, i);
521 return (lb != j+n && *lb == i)
522 ? const_iterator(p,j,lb-j)
523 : end();
524 }
525
526 //===== sizes
527
529 size_type size () const
530 {
531 return n;
532 }
533
534 protected:
536 compressed_base_array_unmanaged ()
537 : n(0), p(0), j(0)
538 {}
539
540 size_type n; // number of elements in array
541 B *p; // pointer to dynamically allocated built-in array
542 size_type* j; // the index set
543 };
544
545} // end namespace Imp
546
547} // end namespace
548
549#endif
#define DUNE_THROW(E, m)
Definition: exceptions.hh:218
constexpr auto equals(T1 &&t1, T2 &&t2)
Equality comparison.
Definition: hybridutilities.hh:587
constexpr auto min
Function object that returns the smaller of the given values.
Definition: hybridutilities.hh:506
constexpr decltype(auto) elementAt(Container &&c, Index &&i)
Get element at given position from container.
Definition: hybridutilities.hh:126
This file implements iterator facade classes for writing stl conformant iterators.
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.111.3 (Nov 21, 23:30, 2024)