Reference documentation for deal.II version 9.1.1
|
#include <deal.II/base/index_set.h>
Classes | |
class | ElementIterator |
class | IntervalAccessor |
class | IntervalIterator |
struct | Range |
Public Types | |
using | size_type = types::global_dof_index |
using | value_type = signed int |
Public Member Functions | |
IndexSet () | |
IndexSet (const size_type size) | |
IndexSet (const IndexSet &)=default | |
IndexSet & | operator= (const IndexSet &)=default |
IndexSet (IndexSet &&is) noexcept | |
IndexSet & | operator= (IndexSet &&is) noexcept |
IndexSet (const Epetra_BlockMap &map) | |
void | clear () |
void | set_size (const size_type size) |
size_type | size () const |
void | add_range (const size_type begin, const size_type end) |
void | add_index (const size_type index) |
template<typename ForwardIterator > | |
void | add_indices (const ForwardIterator &begin, const ForwardIterator &end) |
void | add_indices (const IndexSet &other, const unsigned int offset=0) |
bool | is_element (const size_type index) const |
bool | is_contiguous () const |
bool | is_empty () const |
bool | is_ascending_and_one_to_one (const MPI_Comm &communicator) const |
size_type | n_elements () const |
size_type | nth_index_in_set (const unsigned int local_index) const |
size_type | index_within_set (const size_type global_index) const |
unsigned int | n_intervals () const |
unsigned int | largest_range_starting_index () const |
void | compress () const |
bool | operator== (const IndexSet &is) const |
bool | operator!= (const IndexSet &is) const |
IndexSet | operator & (const IndexSet &is) const |
IndexSet | get_view (const size_type begin, const size_type end) const |
void | subtract_set (const IndexSet &other) |
size_type | pop_back () |
size_type | pop_front () |
void | fill_index_vector (std::vector< size_type > &indices) const |
template<typename VectorType > | |
void | fill_binary_vector (VectorType &vector) const |
template<class StreamType > | |
void | print (StreamType &out) const |
void | write (std::ostream &out) const |
void | read (std::istream &in) |
void | block_write (std::ostream &out) const |
void | block_read (std::istream &in) |
Epetra_Map | make_trilinos_map (const MPI_Comm &communicator=MPI_COMM_WORLD, const bool overlapping=false) const |
std::size_t | memory_consumption () const |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
Iterators | |
ElementIterator | begin () const |
ElementIterator | at (const size_type global_index) const |
ElementIterator | end () const |
IntervalIterator | begin_intervals () const |
IntervalIterator | end_intervals () const |
Static Public Member Functions | |
static ::ExceptionBase & | ExcIndexNotPresent (size_type arg1) |
Private Member Functions | |
void | do_compress () const |
Private Attributes | |
std::vector< Range > | ranges |
bool | is_compressed |
size_type | index_space_size |
size_type | largest_range |
Threads::Mutex | compress_mutex |
Related Functions | |
(Note that these are not member functions.) | |
IndexSet | complete_index_set (const unsigned int N) |
A class that represents a subset of indices among a larger set. For example, it can be used to denote the set of degrees of freedom within the range \([0,\text{dof\_handler.n\_dofs})\) that belongs to a particular subdomain, or those among all degrees of freedom that are stored on a particular processor in a distributed parallel computation.
This class can represent a collection of half-open ranges of indices as well as individual elements. For practical purposes it also stores the overall range these indices can assume. In other words, you need to specify the size of the index space \([0,\text{size})\) of which objects of this class are a subset.
There are two ways to iterate over the IndexSets: First, begin() and end() allow iteration over individual indices in the set. Second, begin_interval() and end_interval() allow iteration over the half-open ranges as described above.
The data structures used in this class along with a rationale can be found in the Distributed Computing paper.
Definition at line 74 of file index_set.h.
size_type
is the type used for storing the size and the individual entries in the IndexSet.
Definition at line 85 of file index_set.h.
using IndexSet::value_type = signed int |
One can see an IndexSet as a container of size size(), where the elements of the containers are bool values that are either false or true, depending on whether a particular index is an element of the IndexSet or not. In other words, an IndexSet is a bit like a vector in which the elements we store are booleans. In this view, the correct local alias indicating the type of the elements of the vector would then be bool
.
On the other hand, bool
has the disadvantage that it is not a numerical type that, for example, allows multiplication with a double
. In other words, one can not easily use a vector of booleans in a place where other vectors are allowed. Consequently, we declare the type of the elements of such a vector as a signed integer. This uses the fact that in the C++ language, booleans are implicitly convertible to integers. In other words, declaring the type of the elements of the vector as a signed integer is only a small lie, but it is a useful one.
Definition at line 104 of file index_set.h.
|
inline |
Default constructor.
Definition at line 1430 of file index_set.h.
|
inlineexplicit |
Constructor that also sets the overall size of the index range.
Definition at line 1438 of file index_set.h.
|
default |
Copy constructor.
|
inlinenoexcept |
Move constructor. Create a new IndexSet by transferring the internal data of the input set.
Definition at line 1446 of file index_set.h.
|
explicit |
Constructor from a Trilinos Epetra_BlockMap.
Definition at line 68 of file index_set.cc.
Move assignment operator. Transfer the internal data of the input set into the current one.
Definition at line 1463 of file index_set.h.
|
inline |
Remove all indices from this index set. The index set retains its size, however.
Definition at line 1576 of file index_set.h.
|
inline |
Set the maximal size of the indices upon which this object operates.
This function can only be called if the index set does not yet contain any elements. This can be achieved by calling clear(), for example.
Definition at line 1588 of file index_set.h.
|
inline |
Return the size of the index space of which this index set is a subset of.
Note that the result is not equal to the number of indices within this set. The latter information is returned by n_elements().
Definition at line 1600 of file index_set.h.
Add the half-open range \([\text{begin},\text{end})\) to the set of indices represented by this class.
[in] | begin | The first element of the range to be added. |
[in] | end | The past-the-end element of the range to be added. |
Definition at line 98 of file index_set.cc.
|
inline |
Add an individual index to the set of indices.
Definition at line 1619 of file index_set.h.
|
inline |
Add a whole set of indices described by dereferencing every element of the iterator range [begin,end)
.
[in] | begin | Iterator to the first element of range of indices to be added |
[in] | end | The past-the-end iterator for the range of elements to be added. |
begin<=end
needs to be satisfied. Definition at line 1641 of file index_set.h.
void IndexSet::add_indices | ( | const IndexSet & | other, |
const unsigned int | offset = 0 |
||
) |
Add the given IndexSet other
to the current one, constructing the union of *this and other
.
If the offset
argument is nonzero, then every index in other
is shifted by offset
before being added to the current index set. This allows to construct, for example, one index set from several others that are supposed to represent index sets corresponding to different ranges (e.g., when constructing the set of nonzero entries of a block vector from the sets of nonzero elements of the individual blocks of a vector).
This function will generate an exception if any of the (possibly shifted) indices of the other
index set lie outside the range [0,size())
represented by the current object.
Definition at line 379 of file index_set.cc.
|
inline |
Return whether the specified index is an element of the index set.
Definition at line 1665 of file index_set.h.
|
inline |
Return whether the index set stored by this object defines a contiguous range. This is true also if no indices are stored at all.
Definition at line 1715 of file index_set.h.
|
inline |
Return whether the index set stored by this object contains no elements. This is similar, but faster than checking n_elements() == 0
.
Definition at line 1724 of file index_set.h.
bool IndexSet::is_ascending_and_one_to_one | ( | const MPI_Comm & | communicator | ) | const |
Return whether the IndexSets are ascending with respect to MPI process number and 1:1, i.e., each index is contained in exactly one IndexSet (among those stored on the different processes), each process stores contiguous subset of indices, and the index set on process \(p+1\) starts at the index one larger than the last one stored on process \(p\). In case there is only one MPI process, this just means that the IndexSet is complete.
Definition at line 664 of file index_set.cc.
|
inline |
Return the number of elements stored in this index set.
Definition at line 1732 of file index_set.h.
|
inline |
Return the global index of the local index with number local_index
stored in this index set. local_index
obviously needs to be less than n_elements().
Definition at line 1780 of file index_set.h.
|
inline |
Return the how-manyth element of this set (counted in ascending order) global_index
is. global_index
needs to be less than the size(). This function returns numbers::invalid_dof_index if the index global_index
is not actually a member of this index set, i.e. if is_element(global_index) is false.
Definition at line 1821 of file index_set.h.
|
inline |
Each index set can be represented as the union of a number of contiguous intervals of indices, where if necessary intervals may only consist of individual elements to represent isolated members of the index set.
This function returns the minimal number of such intervals that are needed to represent the index set under consideration.
Definition at line 1757 of file index_set.h.
|
inline |
This function returns the local index of the beginning of the largest range.
Definition at line 1766 of file index_set.h.
|
inline |
Compress the internal representation by merging individual elements with contiguous ranges, etc. This function does not have any external effect.
Definition at line 1608 of file index_set.h.
|
inline |
Comparison for equality of index sets. This operation is only allowed if the size of the two sets is the same (though of course they do not have to have the same number of indices).
Definition at line 1869 of file index_set.h.
|
inline |
Comparison for inequality of index sets. This operation is only allowed if the size of the two sets is the same (though of course they do not have to have the same number of indices).
Definition at line 1882 of file index_set.h.
Return the intersection of the current index set and the argument given, i.e. a set of indices that are elements of both index sets. The two index sets must have the same size (though of course they do not have to have the same number of indices).
This command takes an interval [begin, end)
and returns the intersection of the current index set with the interval, shifted to the range [0, end-begin)
.
In other words, the result of this operation is the intersection of the set represented by the current object and the interval [begin, end)
, as seen within the interval [begin, end)
by shifting the result of the intersection operation to the left by begin
. This corresponds to the notion of a view: The interval [begin, end)
is a window through which we see the set represented by the current object.
Definition at line 238 of file index_set.cc.
void IndexSet::subtract_set | ( | const IndexSet & | other | ) |
Remove all elements contained in other
from this set. In other words, if \(x\) is the current object and \(o\) the argument, then we compute \(x \leftarrow x \backslash o\).
Definition at line 267 of file index_set.cc.
IndexSet::size_type IndexSet::pop_back | ( | ) |
Remove and return the last element of the last range. This function throws an exception if the IndexSet is empty.
Definition at line 339 of file index_set.cc.
IndexSet::size_type IndexSet::pop_front | ( | ) |
Remove and return the first element of the first range. This function throws an exception if the IndexSet is empty.
Definition at line 357 of file index_set.cc.
void IndexSet::fill_index_vector | ( | std::vector< size_type > & | indices | ) | const |
Fill the given vector with all indices contained in this IndexSet.
Definition at line 505 of file index_set.cc.
void IndexSet::fill_binary_vector | ( | VectorType & | vector | ) | const |
Fill the given vector with either zero or one elements, providing a binary representation of this index set. The given vector is assumed to already have the correct size.
The given argument is filled with integer values zero and one, using vector.operator[]
. Thus, any object that has such an operator can be used as long as it allows conversion of integers zero and one to elements of the vector. Specifically, this is the case for classes Vector, BlockVector, but also std::vector<bool>, std::vector<int>, and std::vector<double>.
|
inline |
Output a text representation of this IndexSet to the given stream. Used for testing.
Definition at line 1915 of file index_set.h.
void IndexSet::write | ( | std::ostream & | out | ) | const |
Write the IndexSet into a text based file format, that can be read in again using the read() function.
Definition at line 434 of file index_set.cc.
void IndexSet::read | ( | std::istream & | in | ) |
Construct the IndexSet from a text based representation given by the stream in
written by the write() function.
Definition at line 449 of file index_set.cc.
void IndexSet::block_write | ( | std::ostream & | out | ) | const |
Write the IndexSet into a binary, compact representation, that can be read in again using the block_read() function.
Definition at line 471 of file index_set.cc.
void IndexSet::block_read | ( | std::istream & | in | ) |
Construct the IndexSet from a binary representation given by the stream in
written by the write_block() function.
Definition at line 485 of file index_set.cc.
Epetra_Map IndexSet::make_trilinos_map | ( | const MPI_Comm & | communicator = MPI_COMM_WORLD , |
const bool | overlapping = false |
||
) | const |
Given an MPI communicator, create a Trilinos map object that represents a distribution of vector elements or matrix rows in which we will locally store those elements or rows for which we store the index in the current index set, and all the other elements/rows elsewhere on one of the other MPI processes.
The last argument only plays a role if the communicator is a parallel one, distributing computations across multiple processors. In that case, if the last argument is false, then it is assumed that the index sets this function is called with on all processors are mutually exclusive but together enumerate each index exactly once. In other words, if you call this function on two processors, then the index sets this function is called with must together have all possible indices from zero to size()-1, and no index must appear in both index sets. This corresponds, for example, to the case where we want to split the elements of vectors into unique subsets to be stored on different processors – no element should be owned by more than one processor, but each element must be owned by one.
On the other hand, if the second argument is true, then the index sets can be overlapping, and they also do not need to span the whole index set. This is a useful operation if we want to create vectors that not only contain the locally owned indices, but for example also the elements that correspond to degrees of freedom located on ghost cells. Another application of this method is to select a subset of the elements of a vector, e.g. for extracting only certain solution components.
Definition at line 594 of file index_set.cc.
std::size_t IndexSet::memory_consumption | ( | ) | const |
Determine an estimate for the memory consumption (in bytes) of this object.
Definition at line 727 of file index_set.cc.
|
inline |
Write or read the data of this object to or from a stream for the purpose of serialization
Definition at line 1946 of file index_set.h.
|
inline |
Return an iterator that points at the first index that is contained in this IndexSet.
Definition at line 1483 of file index_set.h.
|
inline |
Return an element iterator pointing to the element with global index global_index
or the next larger element if the index is not in the set. This is equivalent to
If there is no element in this IndexSet at or behind global_index
, this method will return end().
Definition at line 1495 of file index_set.h.
|
inline |
Return an iterator that points one after the last index that is contained in this IndexSet.
Definition at line 1546 of file index_set.h.
|
inline |
Return an Iterator that points at the first interval of this IndexSet.
Definition at line 1555 of file index_set.h.
|
inline |
Return an Iterator that points one after the last interval of this IndexSet.
Definition at line 1567 of file index_set.h.
|
private |
Actually perform the compress() operation.
Definition at line 127 of file index_set.cc.
|
related |
Create and return an index set of size \(N\) that contains every single index within this range. In essence, this function returns an index set created by
This function exists so that one can create and initialize index sets that are complete in one step, or so one can write code like
Definition at line 978 of file index_set.h.
|
mutableprivate |
A set of contiguous ranges of indices that make up (part of) this index set. This variable is always kept sorted.
The variable is marked "mutable" so that it can be changed by compress(), though this of course doesn't change anything about the external representation of this index set.
Definition at line 917 of file index_set.h.
|
mutableprivate |
True if compress() has been called after the last change in the set of indices.
The variable is marked "mutable" so that it can be changed by compress(), though this of course doesn't change anything about the external representation of this index set.
Definition at line 927 of file index_set.h.
|
private |
The overall size of the index range. Elements of this index set have to have a smaller number than this value.
Definition at line 933 of file index_set.h.
|
mutableprivate |
This integer caches the index of the largest range in ranges
. This gives O(1)
access to the range with most elements, while general access costs O(log(n_ranges))
. The largest range is needed for the methods is_element()
, index_within_set()
, nth_index_in_set
. In many applications, the largest range contains most elements (the locally owned range), whereas there are only a few other elements (ghosts).
Definition at line 944 of file index_set.h.
|
mutableprivate |
A mutex that is used to synchronize operations of the do_compress() function that is called from many 'const' functions via compress().
Definition at line 950 of file index_set.h.