File: indexset.hh

package info (click to toggle)
dune-common 2.11.0-1
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 6,048 kB
  • sloc: cpp: 54,403; python: 4,136; sh: 1,657; makefile: 17
file content (1159 lines) | stat: -rw-r--r-- 31,895 bytes parent folder | download
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
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105
1106
1107
1108
1109
1110
1111
1112
1113
1114
1115
1116
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129
1130
1131
1132
1133
1134
1135
1136
1137
1138
1139
1140
1141
1142
1143
1144
1145
1146
1147
1148
1149
1150
1151
1152
1153
1154
1155
1156
1157
1158
1159
// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
// vi: set et ts=4 sw=2 sts=2:
// SPDX-FileCopyrightInfo: Copyright © DUNE Project contributors, see file LICENSE.md in module root
// SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception
#ifndef DUNE_COMMON_PARALLEL_INDEXSET_HH
#define DUNE_COMMON_PARALLEL_INDEXSET_HH

#include <algorithm>
#include <cstdint> // for uint32_t
#include <iostream>

#include <dune/common/arraylist.hh>
#include <dune/common/exceptions.hh>
#include <dune/common/parallel/localindex.hh>
#include <dune/common/parallel/mpitraits.hh>

namespace Dune
{
  /** @addtogroup Common_Parallel
   *
   * @{
   */
  /**
   * @file
   * @brief Provides a map between global and local indices.
   * @author Markus Blatt
   */
  // forward declarations

  template<class TG, class TL>
  class IndexPair;

  /**
   * @brief Print an index pair.
   * @param os The outputstream to print to.
   * @param pair The index pair to print.
   */
  template<class TG, class TL>
  std::ostream& operator<<(std::ostream& os, const IndexPair<TG,TL>& pair);

  template<class TG, class TL>
  bool operator==(const IndexPair<TG,TL>&, const IndexPair<TG,TL>&);

  template<class TG, class TL>
  bool operator!=(const IndexPair<TG,TL>&, const IndexPair<TG,TL>&);

  template<class TG, class TL>
  bool operator<(const IndexPair<TG,TL>&, const IndexPair<TG,TL>&);

  template<class TG, class TL>
  bool operator>(const IndexPair<TG,TL>&, const IndexPair<TG,TL>&);

  template<class TG, class TL>
  bool operator<=(const IndexPair<TG,TL>&, const IndexPair<TG,TL>&);

  template<class TG, class TL>
  bool operator >=(const IndexPair<TG,TL>&, const IndexPair<TG,TL>&);

  template<class TG, class TL>
  bool operator==(const IndexPair<TG,TL>&, const TG&);

  template<class TG, class TL>
  bool operator!=(const IndexPair<TG,TL>&, const TG&);

  template<class TG, class TL>
  bool operator<(const IndexPair<TG,TL>&, const TG&);

  template<class TG, class TL>
  bool operator>(const IndexPair<TG,TL>&, const TG&);

  template<class TG, class TL>
  bool operator<=(const IndexPair<TG,TL>&, const TG&);

  template<class TG, class TL>
  bool operator >=(const IndexPair<TG,TL>&, const TG&);

  template<typename T>
  struct MPITraits;

  /**
   * @brief A pair consisting of a global and local index.
   */
  template<class TG, class TL>
  class IndexPair
  {
    friend std::ostream& operator<< <>(std::ostream&, const IndexPair<TG,TL>&);
    friend bool operator==<>(const IndexPair<TG,TL>&, const IndexPair<TG,TL>&);
    friend bool operator!=<>(const IndexPair<TG,TL>&, const IndexPair<TG,TL>&);
    friend bool operator< <>(const IndexPair<TG,TL>&, const IndexPair<TG,TL>&);
    friend bool operator><>(const IndexPair<TG,TL>&, const IndexPair<TG,TL>&);
    friend bool operator<=<>(const IndexPair<TG,TL>&, const IndexPair<TG,TL>&);
    friend bool operator>=<>(const IndexPair<TG,TL>&, const IndexPair<TG,TL>&);
    friend bool operator==<>(const IndexPair<TG,TL>&, const TG &);
    friend bool operator!=<>(const IndexPair<TG,TL>&, const TG &);
    friend bool operator< <>(const IndexPair<TG,TL>&, const TG &);
    friend bool operator> <>(const IndexPair<TG,TL>&, const TG &);
    friend bool operator<=<>(const IndexPair<TG,TL>&, const TG &);
    friend bool operator>=<>(const IndexPair<TG,TL>&, const TG &);
    friend struct MPITraits<IndexPair<TG,TL> >;

  public:
    /**
     * @brief the type of the global index.
     *
     * This type has to provide at least a operator&lt; for sorting.
     */
    typedef TG GlobalIndex;

    /**
     * @brief the type of the local index.
     *
     * This class to provide the following functions:
     * \code
     * LocalIndex operator=(int);
     * operator int() const;
     * LocalIndexState state() const;
     * void setState(LocalIndexState);
     * \endcode
     */
    typedef TL LocalIndex;

    /**
     * @brief Constructs a new Pair.
     *
     * @param global The global index.
     * @param local The local index.
     */
    IndexPair(const GlobalIndex& global, const LocalIndex& local);

    /**
     * @brief Construct a new Pair.
     */
    IndexPair();
    /**
     * @brief Constructs a new Pair.
     *
     * The local index will be 0.
     * @param global The global index.
     */
    IndexPair(const GlobalIndex& global);

    /**
     * @brief Get the global index.
     *
     * @return The global index.
     */
    inline const GlobalIndex& global() const;

    /**
     * @brief Get the local index.
     *
     * @return The local index.
     */
    inline LocalIndex& local();

    /**
     * @brief Get the local index.
     *
     * @return The local index.
     */
    inline const LocalIndex& local() const;

    /**
     * @brief Set the local index.
     *
     * @param index The index to set it to.
     */
    inline void setLocal(int index);
  private:
    /** @brief The global index. */
    GlobalIndex global_;
    /** @brief The local index. */
    LocalIndex local_;
  };

  /**
   * @brief The states the index set can be in.
   * @see ParallelIndexSet::state_
   */
  enum ParallelIndexSetState
  {
    /**
     * @brief The default mode.
     * Indicates that the index set is ready to be used.
     */
    GROUND,
    /**
     * @brief Indicates that the index set is currently being resized.
     */
    RESIZE
    /**
     * @brief Indicates that all previously deleted indices are now deleted.
     *
       CLEAN,
     **
     * @brief Indicates that the index set is currently being reordered.
     *
       REORDER
     */
  };

  /**
   * @brief Exception indicating that the index set is not in the expected state.
   */
  class InvalidIndexSetState : public InvalidStateException {};

  // Forward declaration
  template<class I> class GlobalLookupIndexSet;

  /**
   * @brief Manager class for the mapping between local indices and globally unique indices.
   *
   * The mapping is between a globally unique id and local index. The local index is consecutive
   * and non persistent while the global id might not be consecutive but definitely is persistent.
   */
  template<typename TG, typename TL, int N=100>
  class ParallelIndexSet
  {
    friend class GlobalLookupIndexSet<ParallelIndexSet<TG,TL,N> >;

  public:
    /**
     * @brief the type of the global index.
     * This type has to provide at least a operator&lt; for sorting.
     */
    typedef TG GlobalIndex;

    /**
     * @brief The type of the local index, e.g. ParallelLocalIndex.
     *
     * This class to provide the following functions:
     * \code
     * LocalIndex operator=(int);
     * operator int() const;
     * LocalIndexState state() const;
     * void setState(LocalIndexState);
     * \endcode
     */
    typedef TL LocalIndex;

    /**
     * @brief The type of the pair stored.
     */
    typedef Dune::IndexPair<GlobalIndex,LocalIndex> IndexPair;

    /**
     * @brief The size of the individual arrays in the underlying ArrayList.
     *
     * The default value is 100.
     * @see ArrayList::size
     */
    constexpr static int arraySize = (N>0) ? N : 1;

    /** @brief The iterator over the pairs. */
    class iterator :
      public ArrayList<IndexPair,N>::iterator
    {
      typedef typename ArrayList<IndexPair,N>::iterator
      Father;
      friend class ParallelIndexSet<GlobalIndex,LocalIndex,N>;
    public:
      iterator(ParallelIndexSet<TG,TL,N>& indexSet, const Father& father)
        : Father(father), indexSet_(&indexSet)
      {}

    private:
      /**
       * @brief Mark the index as deleted.
       *
       * The deleted flag will be set in the local index.
       * The index will be removed in the endResize method of the
       * index set.
       *
       * @exception InvalidIndexSetState only when NDEBUG is not defined
       */
      inline void markAsDeleted() const
      {
#ifndef NDEBUG
        if(indexSet_->state_ != RESIZE)
          DUNE_THROW(InvalidIndexSetState, "Indices can only be removed "
                     <<"while in RESIZE state!");
#endif
        Father::operator*().local().setState(DELETED);
      }

      /** @brief The index set we are an iterator of. */
      ParallelIndexSet<TG,TL,N>* indexSet_;

    };



    /** @brief The constant iterator over the pairs. */
    typedef typename
    ArrayList<IndexPair,N>::const_iterator
    const_iterator;

    /**
     * @brief Constructor.
     */
    ParallelIndexSet();

    /**
     * @brief Get the state the index set is in.
     * @return The state of the index set.
     */
    inline const ParallelIndexSetState& state()
    {
      return state_;
    }

    /**
     * @brief Indicate that the index set is to be resized.
     * @exception InvalidState If index set was not in
     * ParallelIndexSetState::GROUND mode.
     */
    void beginResize();

    /**
     * @brief Add an new index to the set.
     *
     * The local index is created by the default constructor.
     * @param global The globally unique id of the index.
     * @exception InvalidState If index set is not in
     * ParallelIndexSetState::RESIZE mode.
     */
    inline void add(const GlobalIndex& global);

    /**
     * @brief Add an new index to the set.
     *
     * @param global The globally unique id of the index.
     * @param local The local index.
     * @exception InvalidState If index set is not in
     * ParallelIndexSetState::RESIZE mode.
     */
    inline void add(const GlobalIndex& global, const LocalIndex& local);

    /**
     * @brief Mark an index as deleted.
     *
     * The index will be deleted during endResize().
     * @param position An iterator at the position we want to delete.
     * @exception InvalidState If index set is not in ParallelIndexSetState::RESIZE mode.
     */
    inline void markAsDeleted(const iterator& position);

    /**
     * @brief Indicate that the resizing finishes.
     *
     * @warning Invalidates all pointers stored to the elements of this index set.
     * The local indices will be ordered
     * according to the global indices:
     * Let \f$(g_i,l_i)_{i=0}^N \f$ be the set of all indices then \f$l_i < l_j\f$
     * if and
     * only if \f$g_i < g_j\f$ for arbitrary \f$i \neq j\f$.
     * @exception InvalidState If index set was not in
     * ParallelIndexSetState::RESIZE mode.
     */
    void endResize();

    /**
     * @brief Find the index pair with a specific global id.
     *
     * This starts a binary search for the entry and therefore has complexity
     * log(N).
     * @param global The globally unique id of the pair.
     * @return The pair of indices for the id.
     * @warning If the global index is not in the set a wrong or even a
     * null reference might be returned. To be save use the throwing alternative at.
     */
    inline IndexPair&
    operator[](const GlobalIndex& global);

    /**
     * @brief Find the index pair with a specific global id.
     *
     * This starts a binary search for the entry and therefore has complexity
     * log(N).
     * @param global The globally unique id of the pair.
     * @return The pair of indices for the id.
     * @exception RangeError Thrown if the global id is not known.
     */
    inline IndexPair&
    at(const GlobalIndex& global);

    /**
     * @brief Find the index pair with a specific global id.
     *
     * This starts a binary search for the entry and therefore has complexity
     * log(N).
     * @param global The globally unique id of the pair.
     * @return The pair of indices for the id.
     * @exception RangeError Thrown if the global id is not known.
     */
    inline bool
    exists (const GlobalIndex& global) const;

    /**
     * @brief Find the index pair with a specific global id.
     *
     * This starts a binary search for the entry and therefore has complexity
     * log(N).
     * @param global The globally unique id of the pair.
     * @return The pair of indices for the id.
     * @warning If the global index is not in the set a wrong or even a
     * null reference might be returned. To be save use the throwing alternative at.
     */
    inline const IndexPair&
    operator[](const GlobalIndex& global) const;

    /**
     * @brief Find the index pair with a specific global id.
     *
     * This starts a binary search for the entry and therefore has complexity
     * log(N).
     * @param global The globally unique id of the pair.
     * @return The pair of indices for the id.
     * @exception RangeError Thrown if the global id is not known.
     */
    inline const IndexPair&
    at(const GlobalIndex& global) const;

    /**
     * @brief Get an iterator over the indices positioned at the first index.
     * @return Iterator over the local indices.
     */
    inline iterator begin();

    /**
     * @brief Get an iterator over the indices positioned after the last index.
     * @return Iterator over the local indices.
     */
    inline iterator end();

    /**
     * @brief Get an iterator over the indices positioned at the first index.
     * @return Iterator over the local indices.
     */
    inline const_iterator begin() const;

    /**
     * @brief Get an iterator over the indices positioned after the last index.
     * @return Iterator over the local indices.
     */
    inline const_iterator end() const;

    /**
     * @brief Renumbers the local index numbers.
     *
     * After this function returns the indices are
     * consecutively numbered beginning from 0. Let
     * $(g_i,l_i)$, $(g_j,l_j)$ be two arbitrary index
     * pairs with $g_i<g_j$ then after renumbering
     * $l_i<l_j$ will hold.
     */
    inline void renumberLocal();

    /**
     * @brief Get the internal sequence number.
     *
     * Is initially 0 is incremented for each resize.
     * @return The sequence number.
     */
    inline int seqNo() const;

    /**
     * @brief Get the total number (public and nonpublic) indices.
     * @return The total number (public and nonpublic) indices.
     */
    inline size_t size() const;

  private:
    /** @brief The index pairs. */
    ArrayList<IndexPair,N> localIndices_;
    /** @brief The new indices for the RESIZE state. */
    ArrayList<IndexPair,N> newIndices_;
    /** @brief The state of the index set. */
    ParallelIndexSetState state_;
    /** @brief Number to keep track of the number of resizes. */
    int seqNo_;
    /** @brief Whether entries were deleted in resize mode. */
    bool deletedEntries_;
    /**
     * @brief Merges the _localIndices and newIndices arrays and creates a new
     * localIndices array.
     */
    inline void merge();
  };


  /**
   * @brief Print an index set.
   * @param os The outputstream to print to.
   * @param indexSet The index set to print.
   */
  template<class TG, class TL, int N>
  std::ostream& operator<<(std::ostream& os, const ParallelIndexSet<TG,TL,N>& indexSet);

  /**
   * @brief Decorates an index set with the possibility to find a global index
   * that is mapped to a specific local.
   *
   */
  template<class I>
  class GlobalLookupIndexSet
  {
  public:
    /**
     * @brief The type of the index set.
     */
    typedef I ParallelIndexSet;

    /**
     * @brief The type of the local index.
     */
    typedef typename ParallelIndexSet::LocalIndex LocalIndex;

    /**
     * @brief The type of the global index.
     */
    typedef typename ParallelIndexSet::GlobalIndex GlobalIndex;

    /**
     * @brief The iterator over the index pairs.
     */
    typedef typename ParallelIndexSet::const_iterator const_iterator;

    typedef Dune::IndexPair<typename I::GlobalIndex, typename I::LocalIndex> IndexPair;

    /**
     * @brief Constructor.
     * @param indexset The index set we want to be able to lookup the corresponding
     * global index of a local index.
     * @param size The number of indices present, i.e. one more than the maximum local index.
     */
    GlobalLookupIndexSet(const ParallelIndexSet& indexset, std::size_t size);

    /**
     * @brief Constructor.
     * @param indexset The index set we want to be able to lookup the corresponding
     * global index of a local index.
     */
    GlobalLookupIndexSet(const ParallelIndexSet& indexset);

    /**
     * @brief Destructor.
     */
    ~GlobalLookupIndexSet();

    /**
     * @brief Find the index pair with a specific global id.
     *
     * This starts a binary search for the entry and therefore has complexity
     * log(N). This method is forwarded to the underlying index set.
     * @param global The globally unique id of the pair.
     * @return The pair of indices for the id.
     * @exception RangeError Thrown if the global id is not known.
     */
    inline const IndexPair&
    operator[](const GlobalIndex& global) const;

    /**
     * @brief Get the index pair corresponding to a local index.
     */
    inline const IndexPair*
    pair(const std::size_t& local) const;

    /**
     * @brief Get an iterator over the indices positioned at the first index.
     * @return Iterator over the local indices.
     */
    inline const_iterator begin() const;

    /**
     * @brief Get an iterator over the indices positioned after the last index.
     * @return Iterator over the local indices.
     */
    inline const_iterator end() const;

    /**
     * @brief Get the internal sequence number.
     *
     * Is initially 0 is incremented for each resize.
     * @return The sequence number.
     */
    inline int seqNo() const;

    /**
     * @brief Get the total number (public and nonpublic) indices.
     * @return The total number (public and nonpublic) indices.
     */
    inline size_t size() const;
  private:
    /**
     * @brief The index set we lookup in.
     */
    const ParallelIndexSet& indexSet_;

    /**
     * @brief The number of indices.
     */
    std::size_t size_;

    /**
     * @brief Array with the positions of the corresponding index pair of the index set.
     */
    std::vector<const IndexPair*> indices_;

  };


  template<typename T>
  struct LocalIndexComparator
  {
    static bool compare([[maybe_unused]] const T& t1, [[maybe_unused]] const T& t2)
    {
      return false;
    }
  };

  template<class TG, class TL>
  struct IndexSetSortFunctor
  {
    bool operator()(const IndexPair<TG,TL>& i1, const IndexPair<TG,TL>& i2)
    {
      return i1.global()<i2.global() || (i1.global()==i2.global() &&
                                         LocalIndexComparator<TL>::compare(i1.local(),
                                                                           i2.local()));
    }
  };



  template<class TG, class TL>
  inline std::ostream& operator<<(std::ostream& os, const IndexPair<TG,TL>& pair)
  {
    os<<"{global="<<pair.global_<<", local="<<pair.local_<<"}";
    return os;
  }

  template<class TG, class TL, int N>
  inline std::ostream& operator<<(std::ostream& os, const ParallelIndexSet<TG,TL,N>& indexSet)
  {
    typedef typename ParallelIndexSet<TG,TL,N>::const_iterator Iterator;
    Iterator end = indexSet.end();
    os<<"{";
    for(Iterator index = indexSet.begin(); index != end; ++index)
      os<<*index<<" ";
    os<<"}";
    return os;

  }

  template<class TG, class TL>
  inline bool operator==(const IndexPair<TG,TL>& a, const IndexPair<TG,TL>& b)
  {
    return a.global_==b.global_;
  }

  template<class TG, class TL>
  inline bool operator!=(const IndexPair<TG,TL>& a, const IndexPair<TG,TL>& b)
  {
    return a.global_!=b.global_;
  }

  template<class TG, class TL>
  inline bool operator<(const IndexPair<TG,TL>& a, const IndexPair<TG,TL>& b)
  {
    return a.global_<b.global_;
  }

  template<class TG, class TL>
  inline bool operator>(const IndexPair<TG,TL>& a, const IndexPair<TG,TL>& b)
  {
    return a.global_>b.global_;
  }

  template<class TG, class TL>
  inline bool operator<=(const IndexPair<TG,TL>& a, const IndexPair<TG,TL>& b)
  {
    return a.global_<=b.global_;
  }

  template<class TG, class TL>
  inline bool operator >=(const IndexPair<TG,TL>& a, const IndexPair<TG,TL>& b)
  {
    return a.global_>=b.global_;
  }

  template<class TG, class TL>
  inline bool operator==(const IndexPair<TG,TL>& a, const TG& b)
  {
    return a.global_==b;
  }

  template<class TG, class TL>
  inline bool operator!=(const IndexPair<TG,TL>& a, const TG& b)
  {
    return a.global_!=b;
  }

  template<class TG, class TL>
  inline bool operator<(const IndexPair<TG,TL>& a, const TG& b)
  {
    return a.global_<b;
  }

  template<class TG, class TL>
  inline bool operator>(const IndexPair<TG,TL>& a, const TG& b)
  {
    return a.global_>b;
  }

  template<class TG, class TL>
  inline bool operator<=(const IndexPair<TG,TL>& a, const TG& b)
  {
    return a.global_<=b;
  }

  template<class TG, class TL>
  inline bool operator >=(const IndexPair<TG,TL>& a, const TG& b)
  {
    return a.global_>=b;
  }

#ifndef DOXYGEN

  template<class TG, class TL>
  IndexPair<TG,TL>::IndexPair(const TG& global, const TL& local)
    : global_(global), local_(local){}

  template<class TG, class TL>
  IndexPair<TG,TL>::IndexPair(const TG& global)
    : global_(global), local_(){}

  template<class TG, class TL>
  IndexPair<TG,TL>::IndexPair()
    : global_(), local_(){}

  template<class TG, class TL>
  inline const TG& IndexPair<TG,TL>::global() const {
    return global_;
  }

  template<class TG, class TL>
  inline TL& IndexPair<TG,TL>::local() {
    return local_;
  }

  template<class TG, class TL>
  inline const TL& IndexPair<TG,TL>::local() const {
    return local_;
  }

  template<class TG, class TL>
  inline void IndexPair<TG,TL>::setLocal(int local){
    local_=local;
  }

  template<class TG, class TL, int N>
  ParallelIndexSet<TG,TL,N>::ParallelIndexSet()
    : state_(GROUND), seqNo_(0), deletedEntries_()
  {}

  template<class TG, class TL, int N>
  void ParallelIndexSet<TG,TL,N>::beginResize()
  {

    // Checks in unproductive code
#ifndef NDEBUG
    if(state_!=GROUND)
      DUNE_THROW(InvalidIndexSetState,
                 "IndexSet has to be in GROUND state, when "
                 << "beginResize() is called!");
#endif

    state_ = RESIZE;
    deletedEntries_ = false;
  }

  template<class TG, class TL, int N>
  inline void ParallelIndexSet<TG,TL,N>::add(const GlobalIndex& global)
  {
    // Checks in unproductive code
#ifndef NDEBUG
    if(state_ != RESIZE)
      DUNE_THROW(InvalidIndexSetState, "Indices can only be added "
                 <<"while in RESIZE state!");
#endif
    newIndices_.push_back(IndexPair(global));
  }

  template<class TG, class TL, int N>
  inline void ParallelIndexSet<TG,TL,N>::add(const TG& global, const TL& local)
  {
    // Checks in unproductive code
#ifndef NDEBUG
    if(state_ != RESIZE)
      DUNE_THROW(InvalidIndexSetState, "Indices can only be added "
                 <<"while in RESIZE state!");
#endif
    newIndices_.push_back(IndexPair(global,local));
  }

  template<class TG, class TL, int N>
  inline void ParallelIndexSet<TG,TL,N>::markAsDeleted(const iterator& global)
  {
    // Checks in unproductive code
#ifndef NDEBUG
    if(state_ != RESIZE)
      DUNE_THROW(InvalidIndexSetState, "Indices can only be removed "
                 <<"while in RESIZE state!");
#endif
    deletedEntries_ = true;

    global.markAsDeleted();
  }

  template<class TG, class TL, int N>
  void ParallelIndexSet<TG,TL,N>::endResize() {
    // Checks in unproductive code
#ifndef NDEBUG
    if(state_ != RESIZE)
      DUNE_THROW(InvalidIndexSetState, "endResize called while not "
                 <<"in RESIZE state!");
#endif

    std::sort(newIndices_.begin(), newIndices_.end(), IndexSetSortFunctor<TG,TL>());
    merge();
    seqNo_++;
    state_ = GROUND;
  }


  template<class TG, class TL, int N>
  inline void ParallelIndexSet<TG,TL,N>::merge(){
    if(localIndices_.size()==0)
    {
      localIndices_=newIndices_;
      newIndices_.clear();
    }
    else if(newIndices_.size()>0 || deletedEntries_)
    {
      ArrayList<IndexPair,N> tempPairs;

      auto old = localIndices_.begin();
      auto added = newIndices_.begin();
      const auto endold = localIndices_.end();
      const auto endadded = newIndices_.end();

      while(old != endold && added!= endadded)
      {
        if(old->local().state()==DELETED) {
          old.eraseToHere();
        }
        else
        {
          if(old->global() < added->global() ||
             (old->global() == added->global()
              && LocalIndexComparator<TL>::compare(old->local(),added->local())))
          {
            tempPairs.push_back(*old);
            old.eraseToHere();
            continue;
          }else
          {
            tempPairs.push_back(*added);
            added.eraseToHere();
          }
        }
      }

      while(old != endold)
      {
        if(old->local().state()!=DELETED) {
          tempPairs.push_back(*old);
        }
        old.eraseToHere();
      }

      while(added!= endadded)
      {
        tempPairs.push_back(*added);
        added.eraseToHere();
      }
      localIndices_ = tempPairs;
    }
  }


  template<class TG, class TL, int N>
  inline const IndexPair<TG,TL>&
  ParallelIndexSet<TG,TL,N>::at(const TG& global) const
  {
    // perform a binary search
    int low=0, high=localIndices_.size()-1, probe=-1;

    while(low<high)
    {
      probe = (high + low) / 2;
      if(global <= localIndices_[probe].global())
        high = probe;
      else
        low = probe+1;
    }

    if(probe==-1)
      DUNE_THROW(RangeError, "No entries!");

    if( localIndices_[low].global() != global)
      DUNE_THROW(RangeError, "Could not find entry of "<<global);
    else
      return localIndices_[low];
  }

  template<class TG, class TL, int N>
  inline const IndexPair<TG,TL>&
  ParallelIndexSet<TG,TL,N>::operator[](const TG& global) const
  {
    // perform a binary search
    int low=0, high=localIndices_.size()-1, probe=-1;

    while(low<high)
    {
      probe = (high + low) / 2;
      if(global <= localIndices_[probe].global())
        high = probe;
      else
        low = probe+1;
    }

    return localIndices_[low];
  }
  template<class TG, class TL, int N>
  inline IndexPair<TG,TL>& ParallelIndexSet<TG,TL,N>::at(const TG& global)
  {
    // perform a binary search
    int low=0, high=localIndices_.size()-1, probe=-1;

    while(low<high)
    {
      probe = (high + low) / 2;
      if(localIndices_[probe].global() >= global)
        high = probe;
      else
        low = probe+1;
    }

    if(probe==-1)
      DUNE_THROW(RangeError, "No entries!");

    if( localIndices_[low].global() != global)
      DUNE_THROW(RangeError, "Could not find entry of "<<global);
    else
      return localIndices_[low];
  }

  template<class TG, class TL, int N>
  inline bool ParallelIndexSet<TG,TL,N>::exists (const TG& global) const
  {
    // perform a binary search
    int low=0, high=localIndices_.size()-1, probe=-1;

    while(low<high)
    {
      probe = (high + low) / 2;
      if(localIndices_[probe].global() >= global)
        high = probe;
      else
        low = probe+1;
    }

    if(probe==-1)
      return false;

    if( localIndices_[low].global() != global)
      return false;
    return true;
  }

  template<class TG, class TL, int N>
  inline IndexPair<TG,TL>& ParallelIndexSet<TG,TL,N>::operator[](const TG& global)
  {
    // perform a binary search
    int low=0, high=localIndices_.size()-1, probe=-1;

    while(low<high)
    {
      probe = (high + low) / 2;
      if(localIndices_[probe].global() >= global)
        high = probe;
      else
        low = probe+1;
    }

    return localIndices_[low];
  }
  template<class TG, class TL, int N>
  inline typename ParallelIndexSet<TG,TL,N>::iterator
  ParallelIndexSet<TG,TL,N>::begin()
  {
    return iterator(*this, localIndices_.begin());
  }


  template<class TG, class TL, int N>
  inline typename ParallelIndexSet<TG,TL,N>::iterator
  ParallelIndexSet<TG,TL,N>::end()
  {
    return iterator(*this,localIndices_.end());
  }

  template<class TG, class TL, int N>
  inline typename ParallelIndexSet<TG,TL,N>::const_iterator
  ParallelIndexSet<TG,TL,N>::begin() const
  {
    return localIndices_.begin();
  }


  template<class TG, class TL, int N>
  inline typename ParallelIndexSet<TG,TL,N>::const_iterator
  ParallelIndexSet<TG,TL,N>::end() const
  {
    return localIndices_.end();
  }

  template<class TG, class TL, int N>
  void ParallelIndexSet<TG,TL,N>::renumberLocal(){
#ifndef NDEBUG
    if(state_==RESIZE)
      DUNE_THROW(InvalidIndexSetState, "IndexSet has to be in "
                 <<"GROUND state for renumberLocal()");
#endif

    const auto end_ = end();
    uint32_t index=0;

    for(auto pair=begin(); pair!=end_; index++, ++pair)
      pair->local()=index;
  }

  template<class TG, class TL, int N>
  inline int ParallelIndexSet<TG,TL,N>::seqNo() const
  {
    return seqNo_;
  }

  template<class TG, class TL, int N>
  inline size_t ParallelIndexSet<TG,TL,N>::size() const
  {
    return localIndices_.size();
  }

  template<class I>
  GlobalLookupIndexSet<I>::GlobalLookupIndexSet(const I& indexset,
                                                std::size_t size)
    : indexSet_(indexset), size_(size),
      indices_(size_, static_cast<const IndexPair*>(0))
  {
    const_iterator end_ = indexSet_.end();

    for(const_iterator pair = indexSet_.begin(); pair!=end_; ++pair) {
      assert(pair->local()<size_);
      indices_[pair->local()] = &(*pair);
    }
  }

  template<class I>
  GlobalLookupIndexSet<I>::GlobalLookupIndexSet(const I& indexset)
    : indexSet_(indexset), size_(0)
  {
    const_iterator end_ = indexSet_.end();
    for(const_iterator pair = indexSet_.begin(); pair!=end_; ++pair)
      size_=std::max<std::size_t>(size_,pair->local());

    indices_.resize(++size_,  0);

    for(const_iterator pair = indexSet_.begin(); pair!=end_; ++pair)
      indices_[pair->local()] = &(*pair);
  }

  template<class I>
  GlobalLookupIndexSet<I>::~GlobalLookupIndexSet()
  {}

  template<class I>
  inline const IndexPair<typename I::GlobalIndex, typename I::LocalIndex>*
  GlobalLookupIndexSet<I>::pair(const std::size_t& local) const
  {
    return indices_[local];
  }

  template<class I>
  inline const IndexPair<typename I::GlobalIndex, typename I::LocalIndex>&
  GlobalLookupIndexSet<I>::operator[](const GlobalIndex& global) const
  {
    return indexSet_[global];
  }

  template<class I>
  typename I::const_iterator GlobalLookupIndexSet<I>::begin() const
  {
    return indexSet_.begin();
  }

  template<class I>
  typename I::const_iterator GlobalLookupIndexSet<I>::end() const
  {
    return indexSet_.end();
  }

  template<class I>
  inline size_t GlobalLookupIndexSet<I>::size() const
  {
    return size_;
  }

  template<class I>
  inline int GlobalLookupIndexSet<I>::seqNo() const
  {
    return indexSet_.seqNo();
  }

  template<typename TG, typename TL, int N, typename TG1, typename TL1, int N1>
  bool operator==(const ParallelIndexSet<TG,TL,N>& idxset,
                  const ParallelIndexSet<TG1,TL1,N1>& idxset1)
  {
    if(idxset.size()!=idxset1.size())
      return false;
    typedef typename ParallelIndexSet<TG,TL,N>::const_iterator Iter;
    typedef typename ParallelIndexSet<TG1,TL1,N1>::const_iterator Iter1;
    Iter iter=idxset.begin();
    for(Iter1 iter1=idxset1.begin(); iter1 != idxset1.end(); ++iter, ++iter1) {
      if(iter1->global()!=iter->global())
        return false;
      typedef typename ParallelIndexSet<TG,TL,N>::LocalIndex PI;
      const PI& pi=iter->local(), pi1=iter1->local();

      if(pi!=pi1)
        return false;
    }
    return true;
  }

  template<typename TG, typename TL, int N, typename TG1, typename TL1, int N1>
  bool operator!=(const ParallelIndexSet<TG,TL,N>& idxset,
                  const ParallelIndexSet<TG1,TL1,N1>& idxset1)
  {
    return !(idxset==idxset1);
  }


#endif // DOXYGEN

}

#endif // DUNE_COMMON_PARALLEL_INDEXSET_HH