16 #ifndef dealii__index_set_h 17 #define dealii__index_set_h 19 #include <deal.II/base/config.h> 20 #include <deal.II/base/utilities.h> 21 #include <deal.II/base/exceptions.h> 22 #include <deal.II/base/thread_management.h> 23 #include <boost/serialization/vector.hpp> 29 #ifdef DEAL_II_WITH_TRILINOS 30 # include <Epetra_Map.h> 33 #if defined(DEAL_II_WITH_MPI) || defined(DEAL_II_WITH_PETSC) 37 # ifndef MPI_COMM_WORLD 38 # define MPI_COMM_WORLD 0 42 DEAL_II_NAMESPACE_OPEN
111 #ifdef DEAL_II_WITH_CXX11 135 #ifdef DEAL_II_WITH_TRILINOS 139 explicit IndexSet(
const Epetra_Map &map);
188 template <
typename ForwardIterator>
190 const ForwardIterator &
end);
208 const unsigned int offset = 0);
354 template <
typename VectorType>
361 template <
class StreamType>
362 void print(StreamType &out)
const;
368 void write(std::ostream &out)
const;
374 void read(std::istream &in);
388 #ifdef DEAL_II_WITH_TRILINOS 418 const bool overlapping =
false)
const;
429 <<
"The global index " << arg1
430 <<
" is not an element of this set.");
436 template <
class Archive>
437 void serialize (Archive &ar,
const unsigned int version);
776 inline bool operator< (
const Range &range_1,
777 const Range &range_2)
779 return ((range_1.begin < range_2.begin)
781 ((range_1.begin == range_2.begin)
783 (range_1.end < range_2.end)));
788 return x.end < y.end;
794 return (x.nth_index_in_set+(x.end-x.begin) <
795 y.nth_index_in_set+(y.end-y.begin));
799 inline bool operator== (
const Range &range_1,
800 const Range &range_2)
802 return ((range_1.begin == range_2.begin)
804 (range_1.end == range_2.end));
807 static std::size_t memory_consumption ()
809 return sizeof(
Range);
816 template <
class Archive>
817 void serialize (Archive &ar,
const unsigned int version);
902 : index_set(idxset), range_idx(range_idx)
909 : index_set(idxset), range_idx(
numbers::invalid_dof_index)
916 return index_set->ranges[range_idx].end - index_set->ranges[range_idx].begin;
922 return index_set != NULL && range_idx < index_set->n_intervals();
938 if (range_idx < index_set->
ranges.size()-1)
941 return index_set->end();
950 return index_set->ranges[range_idx].end-1;
955 : index_set (other.index_set), range_idx(other.range_idx)
987 Assert(is_valid(),
ExcMessage(
"Impossible to advance an IndexSet::IntervalIterator that is invalid"));
991 if (range_idx>=index_set->ranges.size())
999 : accessor(idxset, range_idx)
1014 : accessor(other.accessor)
1039 accessor.advance ();
1066 return !(*
this == other);
1084 return static_cast<int>(lhs - rhs);
1086 return -
static_cast<int>(rhs - lhs);
1098 (range_idx < index_set->
ranges.size() && idx<index_set->ranges[range_idx].end)
1101 return range_idx < index_set->ranges.size() && idx<index_set->ranges[range_idx].end;
1106 : index_set(idxset), range_idx(range_idx), idx(index)
1109 ExcMessage(
"Invalid range index for IndexSet::ElementIterator constructor."));
1113 ExcInternalError(
"Invalid index argument for IndexSet::ElementIterator constructor."));
1118 : index_set(idxset), range_idx(
numbers::invalid_dof_index), idx(
numbers::invalid_dof_index)
1125 Assert(is_valid(),
ExcMessage(
"Impossible to dereference an IndexSet::ElementIterator that is invalid"));
1139 Assert(is_valid(),
ExcMessage(
"Impossible to advance an IndexSet::ElementIterator that is invalid"));
1140 if (idx < index_set->
ranges[range_idx].
end)
1143 if (idx == index_set->ranges[range_idx].end)
1146 if (range_idx < index_set->
ranges.size()-1)
1149 idx = index_set->ranges[range_idx].begin;
1180 return !(*
this == other);
1196 if (!(*
this < other))
1197 return -(other-*
this);
1205 std::ptrdiff_t c = index_set->ranges[range_idx].end-idx;
1208 for (
size_type range=range_idx+1; range<index_set->ranges.size() && range<=other.
range_idx; ++range)
1209 c += index_set->ranges[range].end-index_set->ranges[range].begin;
1212 ExcMessage(
"Inconsistent iterator state. Did you invalidate iterators by modifying the IndexSet?"));
1216 c -= index_set->ranges[other.
range_idx].end - other.
idx;
1260 ExcIndexRangeType<size_type> (global_index, 0,
size()));
1267 Range r (global_index, global_index+1);
1270 std::vector<Range>::const_iterator range_begin, range_end;
1271 if (global_index<main_range->
begin)
1273 range_begin =
ranges.begin();
1274 range_end = main_range;
1278 range_begin = main_range;
1279 range_end =
ranges.end();
1284 const std::vector<Range>::const_iterator
1297 if (global_index < p->
begin)
1344 is_compressed (true),
1345 index_space_size (size),
1346 largest_range (
numbers::invalid_unsigned_int)
1351 #ifdef DEAL_II_WITH_CXX11 1356 ranges (
std::move(is.ranges)),
1357 is_compressed (is.is_compressed),
1358 index_space_size (is.index_space_size),
1359 largest_range (is.largest_range)
1362 is.is_compressed =
true;
1363 is.index_space_size = 0;
1373 ranges = std::move (is.ranges);
1379 is.is_compressed =
true;
1380 is.index_space_size = 0;
1409 ExcMessage (
"This function can only be called if the current " 1410 "object does not yet contain any elements."));
1445 const Range new_range(index, index+1);
1447 ranges.push_back(new_range);
1448 else if (index ==
ranges.back().end)
1460 template <
typename ForwardIterator>
1464 const ForwardIterator &end)
1468 for (ForwardIterator p=
begin; p!=
end;)
1472 ForwardIterator q = p;
1474 while ((q !=
end) && (*q == end_index))
1491 if (
ranges.empty() ==
false)
1513 std::vector<Range>::const_iterator
1522 return ((index >= p->begin) && (index < p->end));
1531 return (p->end > index);
1545 return (
ranges.size() <= 1);
1569 v = r.nth_index_in_set + r.end - r.begin;
1574 for (std::vector<Range>::iterator range =
ranges.begin();
1577 s += (range->end - range->begin);
1605 return main_range->nth_index_in_set;
1621 if (n>=main_range->nth_index_in_set &&
1623 return main_range->begin + (n-main_range->nth_index_in_set);
1629 r.nth_index_in_set = n;
1630 std::vector<Range>::const_iterator range_begin, range_end;
1633 range_begin =
ranges.begin();
1634 range_end = main_range;
1638 range_begin = main_range + 1;
1639 range_end =
ranges.end();
1642 const std::vector<Range>::const_iterator
1644 Range::nth_index_compare);
1647 return p->begin + (n-p->nth_index_in_set);
1658 Assert (n <
size(), ExcIndexRangeType<size_type> (n, 0,
size()));
1668 if (n >= main_range->begin && n < main_range->
end)
1669 return (n-main_range->begin) + main_range->nth_index_in_set;
1672 std::vector<Range>::const_iterator range_begin, range_end;
1673 if (n<main_range->
begin)
1675 range_begin =
ranges.begin();
1676 range_end = main_range;
1680 range_begin = main_range + 1;
1681 range_end =
ranges.end();
1684 std::vector<Range>::const_iterator
1686 Range::end_compare);
1689 if (p==range_end || p->end == n || p->begin > n)
1695 return (n-p->begin) + p->nth_index_in_set;
1730 template <
typename Vector>
1739 std::fill (vector.
begin(), vector.
end(), 0);
1743 for (std::vector<Range>::iterator it =
ranges.begin();
1752 template <
class StreamType>
1759 std::vector<Range>::const_iterator p;
1762 if (p->end-p->begin==1)
1765 out <<
"[" << p->begin <<
"," << p->end-1 <<
"]";
1770 out <<
"}" << std::endl;
1775 template <
class Archive>
1780 ar &begin &end &nth_index_in_set;
1785 template <
class Archive>
1793 DEAL_II_NAMESPACE_CLOSE
Iterator lower_bound(Iterator first, Iterator last, const T &val)
void fill_binary_vector(VectorType &vector) const
static const unsigned int invalid_unsigned_int
ElementIterator end() const
void add_index(const size_type index)
IndexSet operator&(const IndexSet &is) const
const IntervalAccessor * operator->() const
const IntervalAccessor & operator*() const
size_type nth_index_in_set(const unsigned int local_index) const
const IndexSet * index_set
IndexSet & operator=(const IndexSet &)=default
void block_read(std::istream &in)
std::ptrdiff_t operator-(const ElementIterator &p) const
void add_indices(const ForwardIterator &begin, const ForwardIterator &end)
unsigned int largest_range_starting_index() const
ElementIterator(const IndexSet *idxset, const size_type range_idx, const size_type index)
bool operator<(const IntervalIterator &) const
std::size_t memory_consumption() const
size_type n_elements() const
size_type index_space_size
bool is_ascending_and_one_to_one(const MPI_Comm &communicator) const
Epetra_Map make_trilinos_map(const MPI_Comm &communicator=MPI_COMM_WORLD, const bool overlapping=false) const
void print(StreamType &out) const
void serialize(Archive &ar, const unsigned int version)
IndexSet complete_index_set(const unsigned int N)
void serialize(Archive &ar, const unsigned int version)
static ::ExceptionBase & ExcMessage(std::string arg1)
void block_write(std::ostream &out) const
IntervalIterator & operator=(const IntervalIterator &other)
std::vector< Range > ranges
static ::ExceptionBase & ExcIndexNotPresent(size_type arg1)
#define DeclException1(Exception1, type1, outsequence)
unsigned int global_dof_index
void subtract_set(const IndexSet &other)
#define Assert(cond, exc)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
bool operator<(const ElementIterator &) const
size_type index_within_set(const size_type global_index) const
bool operator==(const IndexSet &is) const
IntervalAccessor accessor
IntervalAccessor & operator=(const IntervalAccessor &other)
bool operator==(const IntervalIterator &) const
const IndexSet * index_set
bool operator==(const IntervalAccessor &other) const
unsigned int n_intervals() const
int operator-(const IntervalIterator &p) const
bool operator!=(const IndexSet &is) const
void fill_index_vector(std::vector< size_type > &indices) const
void add_range(const size_type begin, const size_type end)
void set_size(const size_type size)
ElementIterator & operator++()
IntervalAccessor(const IndexSet *idxset, const size_type range_idx)
IndexSet get_view(const size_type begin, const size_type end) const
ElementIterator end() const
IntervalIterator begin_intervals() const
bool is_contiguous() const
Threads::Mutex compress_mutex
IntervalIterator end_intervals() const
bool operator<(const IntervalAccessor &other) const
bool operator!=(const ElementIterator &) const
ElementIterator begin() const
void write(std::ostream &out) const
bool is_element(const size_type index) const
bool operator==(const ElementIterator &) const
ElementIterator at(const size_type global_index) const
const types::global_dof_index invalid_dof_index
types::global_dof_index size_type
ElementIterator begin() const
void read(std::istream &in)
size_type n_elements() const
IntervalIterator & operator++()
static ::ExceptionBase & ExcInternalError()
size_type operator*() const
bool operator!=(const IntervalIterator &) const