File: slice.hh

package info (click to toggle)
dune-istl 2.11.0-1
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 3,504 kB
  • sloc: cpp: 34,844; python: 182; sh: 3; makefile: 3
file content (97 lines) | stat: -rw-r--r-- 4,029 bytes parent folder | download | duplicates (2)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
// SPDX-FileCopyrightText: Copyright © DUNE Project contributors, see file LICENSE.md in module root
// SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception
#ifndef DUNE_PYTHON_ISTL_SLICE_HH
#define DUNE_PYTHON_ISTL_SLICE_HH

#include <cassert>

#include <type_traits>

#include <dune/common/iteratorfacades.hh>

namespace Dune
{

  namespace Python
  {

    // ArraySlice
    // ----------

    template< class A >
    struct ArraySlice
    {
      typedef typename A::member_type member_type;
      typedef typename A::size_type size_type;

      typedef decltype( std::declval< A & >()[ 0 ] ) reference;
      typedef decltype( std::declval< const A & >()[ 0 ] ) const_reference;

    private:
      template< class T >
      struct IteratorImpl
        : public RandomAccessIteratorFacade< IteratorImpl< T >, T >
      {
        friend class IteratorImpl< std::add_const_t< T > >;
        friend class IteratorImpl< std::remove_const_t< T > >;

        friend class RandomAccessIteratorFacade< IteratorImpl< std::add_const_t< T > >, std::add_const_t< T > >;
        friend class RandomAccessIteratorFacade< IteratorImpl< std::remove_const_t< T > >, std::remove_const_t< T > >;

        IteratorImpl () noexcept = default;
        IteratorImpl ( T &array, size_type start, size_type step ) noexcept : array_( &array ), index_( start ), step_( step ) {}
        IteratorImpl ( const IteratorImpl< std::remove_const_t< T > > &other ) noexcept : array_( other.array_ ), index_( other.index_ ), step_( other.step_ ) {}

      private:
        bool equals ( const IteratorImpl &other ) const noexcept { return (index_ == other.index_); }

        std::ptrdiff_t distanceTo ( const IteratorImpl &other ) const noexcept { return (other.index_ - index_) / step_; }

        reference elementAt ( std::ptrdiff_t i ) const noexcept( noexcept( std::declval< T & >()[ 0 ] ) ) { return (*array_)[ index_ ]; }

        void increment () noexcept { index_ += step_; }
        void decrement () noexcept { index_ -= step_; }
        void advance ( std::ptrdiff_t i ) noexcept { index_ += i*step_; }

        T *array_ = nullptr;
        size_type index_ = 0, step_ = 0;
      };

    public:
      typedef IteratorImpl< const A > const_iterator;
      typedef IteratorImpl< A > iterator;

      ArraySlice ( A &array, size_type start, size_type step, size_type size ) noexcept
        : array_( array ), start_( start ), step_( step ), size_( size )
      {}

      reference operator[] ( size_type i ) noexcept( noexcept( std::declval< A & >()[ 0 ] ) ) { return array_[ start_ + i*step_ ]; }
      const_reference operator[] ( size_type i ) const noexcept( noexcept( std::declval< const A & >()[ 0 ] ) ) { return array_[ start_ + i*step_ ]; }

      const_iterator begin () const noexcept { return const_iterator( array_, start_, step_ ); }
      iterator begin () noexcept { return iterator( array_, start_, step_); }

      const_iterator end () const noexcept { return const_iterator( array_, start_ + step_*size_, step_ ); }
      iterator end () noexcept { return iterator( array_, start_ + step_*size_, step_ ); }

      const_iterator beforeEnd () const noexcept { return const_iterator( array_, start_ + step_*(size_-1), step_ ); }
      iterator beforeEnd () noexcept { return iterator( array_, start_ + step_*(size_-1), step_ ); }

      const_iterator beforeBegin () const noexcept { return const_iterator( array_, start_ - step_, step_ ); }
      iterator beforeBegin () noexcept { return iterator( array_, start_ - step_, step_ ); }

      const_iterator find ( size_type i ) const noexcept { return iterator( array_, start_ + step_*std::min( i, size_ ), step_ ); }
      iterator find ( size_type i ) noexcept { return iterator( array_, start_ + step_*std::min( i, size_ ), step_ ); }

      size_type size () const noexcept { return size_; }

    private:
      A &array_;
      size_type start_, step_, size_;
    };

  } // namespace Python

} // namespace Dune

#endif // #ifndef DUNE_PYTHON_ISTL_SLICE_HH