|
Reference documentation for deal.II version 9.2.0
|
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\)
\(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\)
\(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\)
\(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Go to the documentation of this file.
16 #ifndef dealii_index_set_h
17 #define dealii_index_set_h
25 #include <boost/serialization/vector.hpp>
32 #ifdef DEAL_II_WITH_TRILINOS
33 # include <Epetra_Map.h>
34 # ifdef DEAL_II_TRILINOS_WITH_TPETRA
35 # include <Tpetra_Map.hpp>
39 #if defined(DEAL_II_WITH_MPI) || defined(DEAL_II_WITH_PETSC)
43 # ifndef MPI_COMM_WORLD
44 # define MPI_COMM_WORLD 0
141 #ifdef DEAL_II_WITH_TRILINOS
145 explicit IndexSet(
const Epetra_BlockMap &map);
198 template <
typename ForwardIterator>
353 std::vector<IndexSet>
355 const std::vector<types::global_dof_index> &n_indices_per_block)
const;
416 template <
typename VectorType>
424 template <
class StreamType>
426 print(StreamType &out)
const;
433 write(std::ostream &out)
const;
440 read(std::istream &in);
456 #ifdef DEAL_II_WITH_TRILINOS
487 const bool overlapping =
false)
const;
489 # ifdef DEAL_II_TRILINOS_WITH_TPETRA
490 Tpetra::Map<int, types::global_dof_index>
492 const bool overlapping =
false)
const;
506 <<
"The global index " << arg1
507 <<
" is not an element of this set.");
513 template <
class Archive>
515 serialize(Archive &ar,
const unsigned int version);
933 return sizeof(
Range);
940 template <
class Archive>
942 serialize(Archive &ar,
const unsigned int version);
1031 , range_idx(range_idx)
1048 : index_set(other.index_set)
1049 , range_idx(other.range_idx)
1061 return index_set->ranges[range_idx].end - index_set->ranges[range_idx].begin;
1069 return index_set !=
nullptr && range_idx < index_set->n_intervals();
1078 return {index_set, range_idx, index_set->ranges[range_idx].begin};
1089 if (range_idx < index_set->
ranges.size() - 1)
1090 return {index_set, range_idx + 1, index_set->ranges[range_idx + 1].begin};
1092 return index_set->end();
1102 return index_set->ranges[range_idx].end - 1;
1125 "Can not compare accessors pointing to different IndexSets"));
1137 "Can not compare accessors pointing to different IndexSets"));
1149 "Impossible to advance an IndexSet::IntervalIterator that is invalid"));
1153 if (range_idx >= index_set->ranges.size())
1163 : accessor(idxset, range_idx)
1228 return !(*
this == other);
1248 "Can not compare iterators belonging to different IndexSets"));
1251 accessor.index_set->ranges.size() :
1255 accessor.index_set->ranges.size() :
1259 return static_cast<int>(lhs - rhs);
1261 return -
static_cast<int>(rhs - lhs);
1273 , range_idx(range_idx)
1278 "Invalid range index for IndexSet::ElementIterator constructor."));
1283 "Invalid index argument for IndexSet::ElementIterator constructor."));
1301 (range_idx < index_set->
ranges.size() &&
1302 idx < index_set->
ranges[range_idx].end),
1305 return (range_idx < index_set->
ranges.size() &&
1306 idx < index_set->
ranges[range_idx].end);
1316 "Impossible to dereference an IndexSet::ElementIterator that is invalid"));
1328 "Can not compare iterators belonging to different IndexSets"));
1329 return range_idx == other.
range_idx && idx == other.
idx;
1340 "Impossible to advance an IndexSet::ElementIterator that is invalid"));
1341 if (idx < index_set->
ranges[range_idx].
end)
1344 if (idx == index_set->ranges[range_idx].end)
1347 if (range_idx < index_set->
ranges.size() - 1)
1350 idx = index_set->ranges[range_idx].begin;
1386 return !(*
this == other);
1397 "Can not compare iterators belonging to different IndexSets"));
1404 inline std::ptrdiff_t
1410 "Can not compare iterators belonging to different IndexSets"));
1413 if (!(*
this < other))
1414 return -(other - *
this);
1423 std::ptrdiff_t c = index_set->ranges[range_idx].end - idx;
1427 range < index_set->ranges.size() && range <= other.
range_idx;
1429 c += index_set->ranges[range].end - index_set->ranges[range].begin;
1432 other.
range_idx < index_set->ranges.size() ||
1435 "Inconsistent iterator state. Did you invalidate iterators by modifying the IndexSet?"));
1440 c -= index_set->ranges[other.
range_idx].end - other.
idx;
1475 : is_compressed(true)
1476 , index_space_size(size)
1483 : ranges(std::move(is.ranges))
1484 , is_compressed(is.is_compressed)
1485 , index_space_size(is.index_space_size)
1486 , largest_range(is.largest_range)
1489 is.is_compressed =
true;
1490 is.index_space_size = 0;
1501 ranges = std::move(is.ranges);
1502 is_compressed = is.is_compressed;
1503 index_space_size = is.index_space_size;
1504 largest_range = is.largest_range;
1507 is.is_compressed =
true;
1508 is.index_space_size = 0;
1523 return {
this, 0,
ranges[0].begin};
1539 std::vector<Range>::const_iterator main_range =
1545 std::vector<Range>::const_iterator range_begin, range_end;
1546 if (global_index < main_range->
begin)
1548 range_begin =
ranges.begin();
1549 range_end = main_range;
1553 range_begin = main_range;
1554 range_end =
ranges.end();
1559 const std::vector<Range>::const_iterator p =
1572 if (global_index < p->
begin)
1626 ExcMessage(
"This function can only be called if the current "
1627 "object does not yet contain any elements."));
1658 const Range new_range(index, index + 1);
1660 ranges.push_back(new_range);
1661 else if (index ==
ranges.back().end)
1690 ranges.push_back(new_range);
1702 template <
typename ForwardIterator>
1714 std::vector<std::pair<size_type, size_type>> tmp_ranges;
1715 bool ranges_are_sorted =
true;
1716 for (ForwardIterator p =
begin; p !=
end;)
1720 ForwardIterator q = p;
1722 while ((q !=
end) && (*q == end_index))
1728 tmp_ranges.emplace_back(begin_index, end_index);
1735 if (p !=
end && *p < end_index)
1736 ranges_are_sorted =
false;
1739 if (!ranges_are_sorted)
1740 std::sort(tmp_ranges.begin(), tmp_ranges.end());
1750 if (tmp_ranges.size() > 9)
1753 tmp_set.
ranges.reserve(tmp_ranges.size());
1754 for (
const auto &i : tmp_ranges)
1759 for (
const auto &i : tmp_ranges)
1768 if (
ranges.empty() ==
false)
1790 std::vector<Range>::const_iterator p = std::upper_bound(
1798 return ((index >= p->begin) && (index < p->
end));
1806 return (p->end > index);
1819 return (
ranges.size() <= 1);
1847 for (
const auto &range :
ranges)
1848 s += (range.end - range.begin);
1872 const std::vector<Range>::const_iterator main_range =
1875 return main_range->nth_index_in_set;
1889 std::vector<Range>::const_iterator main_range =
1891 if (n >= main_range->nth_index_in_set &&
1893 return main_range->begin + (n - main_range->nth_index_in_set);
1900 std::vector<Range>::const_iterator range_begin, range_end;
1903 range_begin =
ranges.begin();
1904 range_end = main_range;
1908 range_begin = main_range + 1;
1909 range_end =
ranges.end();
1912 const std::vector<Range>::const_iterator p =
1916 return p->begin + (n - p->nth_index_in_set);
1936 std::vector<Range>::const_iterator main_range =
1938 if (n >= main_range->begin && n < main_range->
end)
1939 return (n - main_range->begin) + main_range->nth_index_in_set;
1942 std::vector<Range>::const_iterator range_begin, range_end;
1943 if (n < main_range->
begin)
1945 range_begin =
ranges.begin();
1946 range_end = main_range;
1950 range_begin = main_range + 1;
1951 range_end =
ranges.end();
1954 std::vector<Range>::const_iterator p =
1958 if (p == range_end || p->end == n || p->begin > n)
1964 return (n - p->begin) + p->nth_index_in_set;
1995 template <
typename Vector>
2003 std::fill(vector.begin(), vector.end(), 0);
2007 for (
const auto &range :
ranges)
2008 for (
size_type i = range.begin; i < range.end; ++i)
2014 template <
class StreamType>
2020 std::vector<Range>::const_iterator p;
2023 if (p->end - p->begin == 1)
2026 out <<
"[" << p->begin <<
"," << p->end - 1 <<
"]";
2031 out <<
"}" << std::endl;
2036 template <
class Archive>
2045 template <
class Archive>
std::size_t memory_consumption() const
static std::size_t memory_consumption()
std::ptrdiff_t operator-(const ElementIterator &p) const
void subtract_set(const IndexSet &other)
const types::global_dof_index invalid_dof_index
bool operator==(const IndexSet &is) const
const IndexSet * index_set
bool is_contiguous() const
ElementIterator end() const
void write(std::ostream &out) const
bool operator==(const IntervalAccessor &other) const
friend bool operator==(const Range &range_1, const Range &range_2)
#define AssertIndexRange(index, range)
IntervalAccessor(const IndexSet *idxset, const size_type range_idx)
void block_write(std::ostream &out) const
static bool end_compare(const IndexSet::Range &x, const IndexSet::Range &y)
IndexSet get_view(const size_type begin, const size_type end) const
ElementIterator begin() const
Epetra_Map make_trilinos_map(const MPI_Comm &communicator=MPI_COMM_WORLD, const bool overlapping=false) const
static bool nth_index_compare(const IndexSet::Range &x, const IndexSet::Range &y)
void serialize(Archive &ar, const unsigned int version)
void read(std::istream &in)
bool operator!=(const ElementIterator &) const
IntervalIterator begin_intervals() const
ElementIterator at(const size_type global_index) const
size_type n_elements() const
void set_size(const size_type size)
bool operator<(const IntervalAccessor &other) const
Threads::Mutex compress_mutex
IntervalIterator end_intervals() const
std::forward_iterator_tag iterator_category
std::ptrdiff_t difference_type
IntervalIterator & operator=(const IntervalIterator &other)=default
IndexSet tensor_product(const IndexSet &other) const
ElementIterator(const IndexSet *idxset, const size_type range_idx, const size_type index)
static ::ExceptionBase & ExcMessage(std::string arg1)
unsigned int global_dof_index
VectorType::value_type * begin(VectorType &V)
const IntervalAccessor * operator->() const
Iterator lower_bound(Iterator first, Iterator last, const T &val)
int operator-(const IntervalIterator &p) const
#define DEAL_II_NAMESPACE_OPEN
Tpetra::Map< int, types::global_dof_index > make_tpetra_map(const MPI_Comm &communicator=MPI_COMM_WORLD, const bool overlapping=false) const
bool is_ascending_and_one_to_one(const MPI_Comm &communicator) const
bool operator<(const IntervalIterator &) const
bool is_element(const size_type index) const
void fill_binary_vector(VectorType &vector) const
VectorType::value_type * end(VectorType &V)
#define DeclException1(Exception1, type1, outsequence)
static ::ExceptionBase & ExcIndexNotPresent(size_type arg1)
types::global_dof_index size_type
ElementIterator end() const
size_type largest_range_starting_index() const
ElementIterator begin() const
bool operator==(const IntervalIterator &) const
size_type index_space_size
IntervalAccessor accessor
ElementIterator & operator++()
void advance(std::tuple< I1, I2 > &t, const unsigned int n)
IndexSet operator&(const IndexSet &is) const
static ::ExceptionBase & ExcInternalError()
bool operator<(const ElementIterator &) const
#define Assert(cond, exc)
void serialize(Archive &ar, const unsigned int version)
size_type nth_index_in_set(const size_type local_index) const
bool operator!=(const IntervalIterator &) const
bool operator==(const ElementIterator &) const
IndexSet & operator=(const IndexSet &)=default
void add_index(const size_type index)
size_type operator*() const
void block_read(std::istream &in)
static const unsigned int invalid_unsigned_int
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
std::ptrdiff_t difference_type
size_type n_elements() const
size_type index_within_set(const size_type global_index) const
size_type nth_index_in_set
bool operator!=(const IndexSet &is) const
const IntervalAccessor & operator*() const
void fill_index_vector(std::vector< size_type > &indices) const
#define DEAL_II_NAMESPACE_CLOSE
void add_indices(const ForwardIterator &begin, const ForwardIterator &end)
std::vector< IndexSet > split_by_block(const std::vector< types::global_dof_index > &n_indices_per_block) const
std::string compress(const std::string &input)
TrilinosWrappers::types::int_type global_index(const Epetra_BlockMap &map, const ::types::global_dof_index i)
IntervalIterator & operator++()
void print(StreamType &out) const
unsigned int n_intervals() const
std::forward_iterator_tag iterator_category
const IndexSet * index_set
friend bool operator<(const Range &range_1, const Range &range_2)
IndexSet complete_index_set(const IndexSet::size_type N)
std::vector< Range > ranges
void add_range(const size_type begin, const size_type end)
IntervalAccessor & operator=(const IntervalAccessor &other)