Reference documentation for deal.II version GIT relicensing-136-gb80d0be4af 2024-03-18 08:20:02+00:00
\(\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\}}\)
Loading...
Searching...
No Matches
Classes | Public Types | Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes | Related Symbols | List of all members
IndexSet Class Reference

#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
 
IndexSetoperator= (const IndexSet &)=default
 
 IndexSet (IndexSet &&is) noexcept
 
IndexSetoperator= (IndexSet &&is) noexcept
 
 IndexSet (Teuchos::RCP< const Tpetra::Map< int, types::signed_global_dof_index > > map)
 
 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 size_type 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 size_type local_index) const
 
size_type index_within_set (const size_type global_index) const
 
unsigned int n_intervals () const
 
size_type 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
 
IndexSet get_view (const IndexSet &mask) const
 
std::vector< IndexSetsplit_by_block (const std::vector< types::global_dof_index > &n_indices_per_block) const
 
void subtract_set (const IndexSet &other)
 
IndexSet tensor_product (const IndexSet &other) const
 
size_type pop_back ()
 
size_type pop_front ()
 
std::vector< size_typeget_index_vector () const
 
void fill_index_vector (std::vector< size_type > &indices) const
 
template<typename VectorType >
void fill_binary_vector (VectorType &vector) const
 
bool is_subset_of (const IndexSet &other) const
 
template<typename 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
 
Tpetra::Map< int, types::signed_global_dof_indexmake_tpetra_map (const MPI_Comm communicator=MPI_COMM_WORLD, const bool overlapping=false) const
 
Teuchos::RCP< Tpetra::Map< int, types::signed_global_dof_index > > make_tpetra_map_rcp (const MPI_Comm communicator=MPI_COMM_WORLD, const bool overlapping=false) const
 
IS make_petsc_is (const MPI_Comm communicator=MPI_COMM_WORLD) const
 
std::size_t memory_consumption () const
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 
template<typename Vector >
void fill_binary_vector (Vector &vector) const
 
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 ::ExceptionBaseExcIndexNotPresent (size_type arg1)
 

Private Member Functions

void do_compress () const
 
bool is_element_binary_search (const size_type local_index) const
 
size_type nth_index_in_set_binary_search (const size_type local_index) const
 
size_type index_within_set_binary_search (const size_type global_index) const
 
void add_range_lower_bound (const Range &range)
 
void add_ranges_internal (boost::container::small_vector< std::pair< size_type, size_type >, 200 > &tmp_ranges, const bool ranges_are_sorted)
 

Private Attributes

std::vector< Rangeranges
 
bool is_compressed
 
size_type index_space_size
 
size_type largest_range
 
Threads::Mutex compress_mutex
 

Related Symbols

(Note that these are not member symbols.)

IndexSet complete_index_set (const IndexSet::size_type N)
 

Detailed Description

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,\mathrm{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 66 of file index_set.h.

Member Typedef Documentation

◆ size_type

size_type is the type used for storing the size and the individual entries in the IndexSet.

Definition at line 77 of file index_set.h.

◆ value_type

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 96 of file index_set.h.

Constructor & Destructor Documentation

◆ IndexSet() [1/6]

IndexSet::IndexSet ( )
inline

Default constructor.

Definition at line 1646 of file index_set.h.

◆ IndexSet() [2/6]

IndexSet::IndexSet ( const size_type  size)
inlineexplicit

Constructor that also sets the overall size of the index range.

Definition at line 1654 of file index_set.h.

◆ IndexSet() [3/6]

IndexSet::IndexSet ( const IndexSet )
default

Copy constructor.

◆ IndexSet() [4/6]

IndexSet::IndexSet ( IndexSet &&  is)
inlinenoexcept

Move constructor. Create a new IndexSet by transferring the internal data of the input set.

Definition at line 1662 of file index_set.h.

◆ IndexSet() [5/6]

IndexSet::IndexSet ( Teuchos::RCP< const Tpetra::Map< int, types::signed_global_dof_index > >  map)
explicit

Constructor from a Trilinos Teuchos::RCP<Tpetra::Map>.

Definition at line 42 of file index_set.cc.

◆ IndexSet() [6/6]

IndexSet::IndexSet ( const Epetra_BlockMap &  map)
explicit

Constructor from a Trilinos Epetra_BlockMap.

Definition at line 106 of file index_set.cc.

Member Function Documentation

◆ operator=() [1/2]

IndexSet & IndexSet::operator= ( const IndexSet )
default

Copy assignment operator.

◆ operator=() [2/2]

IndexSet & IndexSet::operator= ( IndexSet &&  is)
inlinenoexcept

Move assignment operator. Transfer the internal data of the input set into the current one.

Definition at line 1679 of file index_set.h.

◆ clear()

void IndexSet::clear ( )
inline

Remove all indices from this index set. The index set retains its size, however.

Definition at line 1741 of file index_set.h.

◆ set_size()

void IndexSet::set_size ( const size_type  size)
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 1753 of file index_set.h.

◆ size()

IndexSet::size_type IndexSet::size ( ) const
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 1765 of file index_set.h.

◆ add_range()

void IndexSet::add_range ( const size_type  begin,
const size_type  end 
)
inline

Add the half-open range \([\text{begin},\text{end})\) to the set of indices represented by this class.

Parameters
[in]beginThe first element of the range to be added.
[in]endThe past-the-end element of the range to be added.

Definition at line 1792 of file index_set.h.

◆ add_index()

void IndexSet::add_index ( const size_type  index)
inline

Add an individual index to the set of indices.

If you have many indices to add to the set, consider calling the add_indices() function below. It is considerably more efficient, particularly if the indices to be added are stored in an array that is already sorted.

Definition at line 1784 of file index_set.h.

◆ add_indices() [1/2]

template<typename ForwardIterator >
void IndexSet::add_indices ( const ForwardIterator &  begin,
const ForwardIterator &  end 
)
inline

Add a whole set of indices described by dereferencing every element of the iterator range [begin,end).

Parameters
[in]beginIterator to the first element of range of indices to be added.
[in]endThe past-the-end iterator for the range of elements to be added.
Precondition
The condition begin<=end needs to be satisfied. They also obviously have to point into the same container.
Note
The operations of this function are substantially more efficient if the indices pointed to by the range of iterators are already sorted. As a consequence, it is often highly beneficial to sort the range of indices pointed to by the iterator range given by the arguments before calling this function. Note also that the function deals efficiently with sorted indices that contain duplicates (e.g., if the iterator range points to the list (1,1,1,2,2,4,6,8,8,8)). In other words, it is very useful to call std::sort() on the range of indices you are about to add, but it is not necessary to call the usual combination of std::unique and std::erase to reduce the list of indices to only a set of unique elements.

Definition at line 1820 of file index_set.h.

◆ add_indices() [2/2]

void IndexSet::add_indices ( const IndexSet other,
const size_type  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 649 of file index_set.cc.

◆ is_element()

bool IndexSet::is_element ( const size_type  index) const
inline

Return whether the specified index is an element of the index set.

Definition at line 1883 of file index_set.h.

◆ is_contiguous()

bool IndexSet::is_contiguous ( ) const
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 1906 of file index_set.h.

◆ is_empty()

bool IndexSet::is_empty ( ) const
inline

Return whether the index set stored by this object contains no elements. This is similar, but faster than checking n_elements() == 0.

Note that a set being empty does not imply that the size of the index space must be zero. Rather, this function returns true if the subset of indices in the index space is the empty set, regardless of the size of the index space.

Definition at line 1915 of file index_set.h.

◆ is_ascending_and_one_to_one()

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 1125 of file index_set.cc.

◆ n_elements()

IndexSet::size_type IndexSet::n_elements ( ) const
inline

Return the number of elements stored in this index set.

Definition at line 1923 of file index_set.h.

◆ nth_index_in_set()

IndexSet::size_type IndexSet::nth_index_in_set ( const size_type  local_index) const
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 1971 of file index_set.h.

◆ index_within_set()

IndexSet::size_type IndexSet::index_within_set ( const size_type  global_index) const
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 1990 of file index_set.h.

◆ n_intervals()

unsigned int IndexSet::n_intervals ( ) const
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 1948 of file index_set.h.

◆ largest_range_starting_index()

IndexSet::size_type IndexSet::largest_range_starting_index ( ) const
inline

This function returns the local index of the beginning of the largest range.

In other words, the return value is nth_index_in_set(x), where x is the first index of the largest contiguous range of indices in the IndexSet. The return value is therefore equal to the number of elements in the set that come before the largest range.

This call assumes that the IndexSet is nonempty.

Definition at line 1957 of file index_set.h.

◆ compress()

void IndexSet::compress ( ) const
inline

Compress the internal representation by merging individual elements with contiguous ranges, etc. This function does not have any external effect.

Definition at line 1773 of file index_set.h.

◆ operator==()

bool IndexSet::operator== ( const IndexSet is) const
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), or if one of the two objects being compared is empty. The comparison between two objects of different sizes would of course, intuitively, result in a false outcome, but it is often a sign of a programming mistake to compare index sets of different sizes against each other, and the comparison is consequently not allowed. On the other hand, the comparison against an empty object makes sense to ensure, for example, that an IndexSet object has been initialized.

Definition at line 2016 of file index_set.h.

◆ operator!=()

bool IndexSet::operator!= ( const IndexSet is) const
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), or if one of the two objects being compared is empty. The comparison between two objects of different sizes would of course, intuitively, result in a false outcome, but it is often a sign of a programming mistake to compare index sets of different sizes against each other, and the comparison is consequently not allowed. On the other hand, the comparison against an empty object makes sense to ensure, for example, that an IndexSet object has been initialized.

Definition at line 2037 of file index_set.h.

◆ operator&()

IndexSet IndexSet::operator& ( const IndexSet is) const

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).

◆ get_view() [1/2]

IndexSet IndexSet::get_view ( const size_type  begin,
const size_type  end 
) const

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.

A more general function of the same name, taking a mask instead of just an interval to define the view, is below.

Definition at line 270 of file index_set.cc.

◆ get_view() [2/2]

IndexSet IndexSet::get_view ( const IndexSet mask) const

This command takes a "mask", i.e., a second index set of same size as the current one and returns the intersection of the current index set the mask, shifted to the index of an entry within the given mask. For example, if the current object is a an IndexSet object representing an index space [0,100) containing indices [20,40), and if the mask represents an index space of the same size but containing all 50 odd indices in this range, then the result will be an index set for a space of size 50 that contains those indices that correspond to the question "the how many'th entry in the mask are the indices [20,40). This will result in an index set of size 50 that contains the indices {11,12,13,14,15,16,17,18,19,20} (because, for example, the index 20 in the original set is not in the mask, but 21 is and corresponds to the 11th entry of the mask – the mask contains the elements {1,3,5,7,9,11,13,15,17,19,21,...}).

In other words, the result of this operation is the intersection of the set represented by the current object and the mask, as seen within the mask. This corresponds to the notion of a view: The mask is a window through which we see the set represented by the current object.

A typical case where this function is useful is as follows. Say, you have a block linear system in which you have blocks corresponding to variables \((u,p,T,c)\) (which you can think of as velocity, pressure, temperature, and chemical composition – or whatever other kind of problem you are currently considering in your own work). We solve this in parallel, so every MPI process has its own locally_owned_dofs index set that describes which among all \(N_\text{dofs}\) degrees of freedom this process owns. Let's assume we have developed a linear solver or preconditioner that first solves the coupled \(u\)- \(T\) system, and once that is done, solves the \(p\)- \(c\) system. In this case, it is often useful to set up block vectors with only two components corresponding to the \(u\) and \(T\) components, and later for only the \(p\)- \(c\) components of the solution. The question is which of the components of these 2-block vectors are locally owned? The answer is that we need to get a view of the locally_owned_dofs index set in which we apply a mask that corresponds to the variables we're currently interested in. For the \(u\)- \(T\) system, we need a mask (corresponding to an index set of size \(N_\text{dofs}\)) that contains all indices of \(u\) degrees of freedom as well as all indices of \(T\) degrees of freedom. The resulting view is an index set of size \(N_u+N_T\) that contains the indices of the locally owned \(u\) and \(T\) degrees of freedom.

Definition at line 305 of file index_set.cc.

◆ split_by_block()

std::vector< IndexSet > IndexSet::split_by_block ( const std::vector< types::global_dof_index > &  n_indices_per_block) const

Split the set indices represented by this object into blocks given by the n_indices_per_block structure. The sum of its entries must match the global size of the current object.

Definition at line 441 of file index_set.cc.

◆ subtract_set()

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 470 of file index_set.cc.

◆ tensor_product()

IndexSet IndexSet::tensor_product ( const IndexSet other) const

Return a new IndexSet, with global size equal to this->size()*other.size(), containing for every element n of this IndexSet, the entries in the half open range [n*other.size(), (n+1)*other.size()) of the other IndexSet.

The name results from the perspective that one starts with an IndexSet and takes the tensor product with another IndexSet with other.size() elements; this results in a matrix of size this->size() times other.size() that has ones in exactly the rows for which this IndexSet contained an index and in the columns for which the other IndexSet contained an index. This matrix is then "unrolled" again by going through each row one by one and reindexing the entries of the matrix in consecutive order. A one in the matrix then corresponds to an entry in the reindexed IndexSet that is returned by this function.

Definition at line 539 of file index_set.cc.

◆ pop_back()

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 551 of file index_set.cc.

◆ pop_front()

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 569 of file index_set.cc.

◆ get_index_vector()

std::vector< IndexSet::size_type > IndexSet::get_index_vector ( ) const

Return a vector with all indices contained in this IndexSet. This vector may of course be quite large if the IndexSet stores many indices. (This may be true even if the IndexSet itself does not take up much memory: IndexSet stores indices in a compressed format in which contiguous ranges of indices are only stored using pairs of indices.)

Definition at line 923 of file index_set.cc.

◆ fill_index_vector()

void IndexSet::fill_index_vector ( std::vector< size_type > &  indices) const

Fill the given vector with all indices contained in this IndexSet.

This function is equivalent to calling get_index_vector() and assigning the result to the indices argument.

Definition at line 942 of file index_set.cc.

◆ fill_binary_vector() [1/2]

template<typename VectorType >
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>.

◆ is_subset_of()

bool IndexSet::is_subset_of ( const IndexSet other) const

Determine whether the current object represents a set of indices that is a subset of the set represented by the argument. This function returns true if the two sets are the same, that is, it considers the "subset" comparison typically used in set theory, rather than the "strict subset" comparison.

Definition at line 703 of file index_set.cc.

◆ print()

template<typename StreamType >
void IndexSet::print ( StreamType &  out) const
inline

Output a text representation of this IndexSet to the given stream. Used for testing.

Definition at line 2078 of file index_set.h.

◆ write()

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 727 of file index_set.cc.

◆ read()

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 742 of file index_set.cc.

◆ block_write()

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 764 of file index_set.cc.

◆ block_read()

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 778 of file index_set.cc.

◆ make_trilinos_map()

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 1036 of file index_set.cc.

◆ make_tpetra_map()

Tpetra::Map< int, types::signed_global_dof_index > IndexSet::make_tpetra_map ( const MPI_Comm  communicator = MPI_COMM_WORLD,
const bool  overlapping = false 
) const

Definition at line 953 of file index_set.cc.

◆ make_tpetra_map_rcp()

Teuchos::RCP< Tpetra::Map< int, types::signed_global_dof_index > > IndexSet::make_tpetra_map_rcp ( const MPI_Comm  communicator = MPI_COMM_WORLD,
const bool  overlapping = false 
) const

Definition at line 962 of file index_set.cc.

◆ make_petsc_is()

IS IndexSet::make_petsc_is ( const MPI_Comm  communicator = MPI_COMM_WORLD) const

Definition at line 1105 of file index_set.cc.

◆ memory_consumption()

std::size_t IndexSet::memory_consumption ( ) const

Determine an estimate for the memory consumption (in bytes) of this object.

Definition at line 1191 of file index_set.cc.

◆ serialize()

template<class Archive >
void IndexSet::serialize ( Archive &  ar,
const unsigned int  version 
)
inline

Write or read the data of this object to or from a stream for the purpose of serialization using the BOOST serialization library.

Definition at line 2109 of file index_set.h.

◆ begin()

IndexSet::ElementIterator IndexSet::begin ( ) const
inline

Return an iterator that points at the first index that is contained in this IndexSet.

Definition at line 1699 of file index_set.h.

◆ at()

IndexSet::ElementIterator IndexSet::at ( const size_type  global_index) const

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

auto p = begin();
while (*p<global_index)
++p;
return p;
ElementIterator begin() const
Definition index_set.h:1699

If there is no element in this IndexSet at or behind global_index, this method will return end().

Definition at line 873 of file index_set.cc.

◆ end()

IndexSet::ElementIterator IndexSet::end ( ) const
inline

Return an iterator that points one after the last index that is contained in this IndexSet.

Definition at line 1711 of file index_set.h.

◆ begin_intervals()

IndexSet::IntervalIterator IndexSet::begin_intervals ( ) const
inline

Return an Iterator that points at the first interval of this IndexSet.

Definition at line 1720 of file index_set.h.

◆ end_intervals()

IndexSet::IntervalIterator IndexSet::end_intervals ( ) const
inline

Return an Iterator that points one after the last interval of this IndexSet.

Definition at line 1732 of file index_set.h.

◆ do_compress()

void IndexSet::do_compress ( ) const
private

Actually perform the compress() operation.

Definition at line 136 of file index_set.cc.

◆ is_element_binary_search()

bool IndexSet::is_element_binary_search ( const size_type  local_index) const
private

Expensive part of is_element() that does a binary search in case we did not find the index in the largest range. Kept separate to avoid pulling in a binary search in the header and make it easy for the compiler to inline the fast path.

Definition at line 798 of file index_set.cc.

◆ nth_index_in_set_binary_search()

IndexSet::size_type IndexSet::nth_index_in_set_binary_search ( const size_type  local_index) const
private

Expensive part of nth_index_in_set() that does the binary search in case we did not find the index in the largest range. Kept separate to avoid using a binary search in the header and make it easy for the compiler to inline the fast path.

Definition at line 834 of file index_set.cc.

◆ index_within_set_binary_search()

IndexSet::size_type IndexSet::index_within_set_binary_search ( const size_type  global_index) const
private

Expensive part of index_within_set() that does the binary search in case we did not find the index in the largest range. Kept separate to avoid using a binary search in the header and make it easy for the compiler to inline the fast path.

Definition at line 851 of file index_set.cc.

◆ add_range_lower_bound()

void IndexSet::add_range_lower_bound ( const Range range)
private

Expensive part of add_index() and add_range(). Defined in separate function to avoid using a binary search in the header and make it easy for the compiler to inline the fast path.

Definition at line 591 of file index_set.cc.

◆ add_ranges_internal()

void IndexSet::add_ranges_internal ( boost::container::small_vector< std::pair< size_type, size_type >, 200 > &  tmp_ranges,
const bool  ranges_are_sorted 
)
private

Expensive part of add_indices().

Definition at line 607 of file index_set.cc.

◆ fill_binary_vector() [2/2]

template<typename Vector >
void IndexSet::fill_binary_vector ( Vector vector) const

Definition at line 2059 of file index_set.h.

Friends And Related Symbol Documentation

◆ complete_index_set()

IndexSet complete_index_set ( const IndexSet::size_type  N)
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

IndexSet is (N);
is.add_range(0, N);

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

if (my_index_set == complete_index_set(my_index_set.size())
...
IndexSet complete_index_set(const IndexSet::size_type N)
Definition index_set.h:1193

Definition at line 1193 of file index_set.h.

Member Data Documentation

◆ ranges

std::vector<Range> IndexSet::ranges
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 1087 of file index_set.h.

◆ is_compressed

bool IndexSet::is_compressed
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 1097 of file index_set.h.

◆ index_space_size

size_type IndexSet::index_space_size
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 1103 of file index_set.h.

◆ largest_range

size_type IndexSet::largest_range
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 1114 of file index_set.h.

◆ compress_mutex

Threads::Mutex IndexSet::compress_mutex
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 1120 of file index_set.h.


The documentation for this class was generated from the following files: