Reference documentation for deal.II version Git c633c6cdfb 2022-01-24 16:30:41 -0500
\(\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\}}\)
index_set.cc
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2005 - 2020 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 #include <deal.II/base/index_set.h>
18 #include <deal.II/base/mpi.h>
19 
20 #include <vector>
21 
22 #ifdef DEAL_II_WITH_TRILINOS
23 # ifdef DEAL_II_WITH_MPI
24 # include <Epetra_MpiComm.h>
25 # endif
26 # include <Epetra_Map.h>
27 # include <Epetra_SerialComm.h>
28 #endif
29 
31 
32 
33 
34 #ifdef DEAL_II_WITH_TRILINOS
35 
36 // the 64-bit path uses a few different names, so put that into a separate
37 // implementation
38 
39 # ifdef DEAL_II_WITH_64BIT_INDICES
40 
41 IndexSet::IndexSet(const Epetra_BlockMap &map)
42  : is_compressed(true)
43  , index_space_size(1 + map.MaxAllGID64())
44  , largest_range(numbers::invalid_unsigned_int)
45 {
46  Assert(map.MinAllGID64() == 0,
47  ExcMessage(
48  "The Epetra_BlockMap does not contain the global index 0, "
49  "which means some entries are not present on any processor."));
50 
51  // For a contiguous map, we do not need to go through the whole data...
52  if (map.LinearMap())
53  add_range(size_type(map.MinMyGID64()), size_type(map.MaxMyGID64() + 1));
54  else
55  {
56  const size_type n_indices = map.NumMyElements();
57  size_type * indices =
58  reinterpret_cast<size_type *>(map.MyGlobalElements64());
59  add_indices(indices, indices + n_indices);
60  }
61  compress();
62 }
63 
64 # else
65 
66 // this is the standard 32-bit implementation
67 
68 IndexSet::IndexSet(const Epetra_BlockMap &map)
69  : is_compressed(true)
70  , index_space_size(1 + map.MaxAllGID())
71  , largest_range(numbers::invalid_unsigned_int)
72 {
73  Assert(map.MinAllGID() == 0,
74  ExcMessage(
75  "The Epetra_BlockMap does not contain the global index 0, "
76  "which means some entries are not present on any processor."));
77 
78  // For a contiguous map, we do not need to go through the whole data...
79  if (map.LinearMap())
80  add_range(size_type(map.MinMyGID()), size_type(map.MaxMyGID() + 1));
81  else
82  {
83  const size_type n_indices = map.NumMyElements();
84  unsigned int * indices =
85  reinterpret_cast<unsigned int *>(map.MyGlobalElements());
86  add_indices(indices, indices + n_indices);
87  }
88  compress();
89 }
90 
91 # endif
92 
93 #endif // ifdef DEAL_II_WITH_TRILINOS
94 
95 
96 
97 void
99 {
100  // we will, in the following, modify mutable variables. this can only
101  // work in multithreaded applications if we lock the data structures
102  // via a mutex, so that users can call 'const' functions from threads
103  // in parallel (and these 'const' functions can then call compress()
104  // which itself calls the current function)
105  std::lock_guard<std::mutex> lock(compress_mutex);
106 
107  // see if any of the contiguous ranges can be merged. do not use
108  // std::vector::erase in-place as it is quadratic in the number of
109  // ranges. since the ranges are sorted by their first index, determining
110  // overlap isn't all that hard
111  std::vector<Range>::iterator store = ranges.begin();
112  for (std::vector<Range>::iterator i = ranges.begin(); i != ranges.end();)
113  {
114  std::vector<Range>::iterator next = i;
115  ++next;
116 
117  size_type first_index = i->begin;
118  size_type last_index = i->end;
119 
120  // see if we can merge any of the following ranges
121  while (next != ranges.end() && (next->begin <= last_index))
122  {
123  last_index = std::max(last_index, next->end);
124  ++next;
125  }
126  i = next;
127 
128  // store the new range in the slot we last occupied
129  *store = Range(first_index, last_index);
130  ++store;
131  }
132  // use a compact array with exactly the right amount of storage
133  if (store != ranges.end())
134  {
135  std::vector<Range> new_ranges(ranges.begin(), store);
136  ranges.swap(new_ranges);
137  }
138 
139  // now compute indices within set and the range with most elements
140  size_type next_index = 0, largest_range_size = 0;
141  for (std::vector<Range>::iterator i = ranges.begin(); i != ranges.end(); ++i)
142  {
143  Assert(i->begin < i->end, ExcInternalError());
144 
145  i->nth_index_in_set = next_index;
146  next_index += (i->end - i->begin);
147  if (i->end - i->begin > largest_range_size)
148  {
149  largest_range_size = i->end - i->begin;
150  largest_range = i - ranges.begin();
151  }
152  }
153  is_compressed = true;
154 
155  // check that next_index is correct. needs to be after the previous
156  // statement because we otherwise will get into an endless loop
157  Assert(next_index == n_elements(), ExcInternalError());
158 }
159 
160 
161 
162 #ifndef DOXYGEN
163 IndexSet
164 IndexSet::operator&(const IndexSet &is) const
165 {
166  Assert(size() == is.size(), ExcDimensionMismatch(size(), is.size()));
167 
168  compress();
169  is.compress();
170 
171  std::vector<Range>::const_iterator r1 = ranges.begin(),
172  r2 = is.ranges.begin();
173  IndexSet result(size());
174 
175  while ((r1 != ranges.end()) && (r2 != is.ranges.end()))
176  {
177  // if r1 and r2 do not overlap at all, then move the pointer that sits
178  // to the left of the other up by one
179  if (r1->end <= r2->begin)
180  ++r1;
181  else if (r2->end <= r1->begin)
182  ++r2;
183  else
184  {
185  // the ranges must overlap somehow
186  Assert(((r1->begin <= r2->begin) && (r1->end > r2->begin)) ||
187  ((r2->begin <= r1->begin) && (r2->end > r1->begin)),
188  ExcInternalError());
189 
190  // add the overlapping range to the result
191  result.add_range(std::max(r1->begin, r2->begin),
192  std::min(r1->end, r2->end));
193 
194  // now move that iterator that ends earlier one up. note that it has
195  // to be this one because a subsequent range may still have a chance
196  // of overlapping with the range that ends later
197  if (r1->end <= r2->end)
198  ++r1;
199  else
200  ++r2;
201  }
202  }
203 
204  result.compress();
205  return result;
206 }
207 #endif
208 
209 
210 
211 IndexSet
213 {
214  Assert(begin <= end,
215  ExcMessage("End index needs to be larger or equal to begin index!"));
216  Assert(end <= size(), ExcMessage("Given range exceeds index set dimension"));
217 
218  IndexSet result(end - begin);
219  std::vector<Range>::const_iterator r1 = ranges.begin();
220 
221  while (r1 != ranges.end())
222  {
223  if ((r1->end > begin) && (r1->begin < end))
224  {
225  result.add_range(std::max(r1->begin, begin) - begin,
226  std::min(r1->end, end) - begin);
227  }
228  else if (r1->begin >= end)
229  break;
230 
231  ++r1;
232  }
233 
234  result.compress();
235  return result;
236 }
237 
238 std::vector<IndexSet>
240  const std::vector<types::global_dof_index> &n_indices_per_block) const
241 {
242  std::vector<IndexSet> partitioned;
243  const unsigned int n_blocks = n_indices_per_block.size();
244 
245  partitioned.reserve(n_blocks);
246  types::global_dof_index start = 0;
247  for (const auto n_block_indices : n_indices_per_block)
248  {
249  partitioned.push_back(this->get_view(start, start + n_block_indices));
250  start += n_block_indices;
251  }
252 
253 #ifdef DEBUG
255  for (const auto &partition : partitioned)
256  {
257  sum += partition.size();
258  }
259  AssertDimension(sum, this->size());
260 #endif
261 
262  return partitioned;
263 }
264 
265 void
267 {
268  compress();
269  other.compress();
270  is_compressed = false;
271 
272 
273  // we save all new ranges to our IndexSet in an temporary vector and
274  // add all of them in one go at the end.
275  std::vector<Range> new_ranges;
276 
277  std::vector<Range>::iterator own_it = ranges.begin();
278  std::vector<Range>::iterator other_it = other.ranges.begin();
279 
280  while (own_it != ranges.end() && other_it != other.ranges.end())
281  {
282  // advance own iterator until we get an overlap
283  if (own_it->end <= other_it->begin)
284  {
285  new_ranges.push_back(*own_it);
286  ++own_it;
287  continue;
288  }
289  // we are done with other_it, so advance
290  if (own_it->begin >= other_it->end)
291  {
292  ++other_it;
293  continue;
294  }
295 
296  // Now own_it and other_it overlap. First save the part of own_it that
297  // is before other_it (if not empty).
298  if (own_it->begin < other_it->begin)
299  {
300  Range r(own_it->begin, other_it->begin);
301  r.nth_index_in_set = 0; // fix warning of unused variable
302  new_ranges.push_back(r);
303  }
304  // change own_it to the sub range behind other_it. Do not delete own_it
305  // in any case. As removal would invalidate iterators, we just shrink
306  // the range to an empty one.
307  own_it->begin = other_it->end;
308  if (own_it->begin > own_it->end)
309  {
310  own_it->begin = own_it->end;
311  ++own_it;
312  }
313 
314  // continue without advancing iterators, the right one will be advanced
315  // next.
316  }
317 
318  // make sure to take over the remaining ranges
319  for (; own_it != ranges.end(); ++own_it)
320  new_ranges.push_back(*own_it);
321 
322  ranges.clear();
323 
324  // done, now add the temporary ranges
325  const std::vector<Range>::iterator end = new_ranges.end();
326  for (std::vector<Range>::iterator it = new_ranges.begin(); it != end; ++it)
327  add_range(it->begin, it->end);
328 
329  compress();
330 }
331 
332 
333 
334 IndexSet
336 {
337  IndexSet set(this->size() * other.size());
338  for (const auto el : *this)
339  set.add_indices(other, el * other.size());
340  set.compress();
341  return set;
342 }
343 
344 
345 
348 {
349  Assert(is_empty() == false,
350  ExcMessage(
351  "pop_back() failed, because this IndexSet contains no entries."));
352 
353  const size_type index = ranges.back().end - 1;
354  --ranges.back().end;
355 
356  if (ranges.back().begin == ranges.back().end)
357  ranges.pop_back();
358 
359  return index;
360 }
361 
362 
363 
366 {
367  Assert(is_empty() == false,
368  ExcMessage(
369  "pop_front() failed, because this IndexSet contains no entries."));
370 
371  const size_type index = ranges.front().begin;
372  ++ranges.front().begin;
373 
374  if (ranges.front().begin == ranges.front().end)
375  ranges.erase(ranges.begin());
376 
377  // We have to set this in any case, because nth_index_in_set is no longer
378  // up to date for all but the first range
379  is_compressed = false;
380 
381  return index;
382 }
383 
384 
385 
386 void
387 IndexSet::add_indices(const IndexSet &other, const size_type offset)
388 {
389  if ((this == &other) && (offset == 0))
390  return;
391 
392  if (other.ranges.size() != 0)
393  {
394  AssertIndexRange(other.ranges.back().end - 1, index_space_size);
395  }
396 
397  compress();
398  other.compress();
399 
400  std::vector<Range>::const_iterator r1 = ranges.begin(),
401  r2 = other.ranges.begin();
402 
403  std::vector<Range> new_ranges;
404  // just get the start and end of the ranges right in this method, everything
405  // else will be done in compress()
406  while (r1 != ranges.end() || r2 != other.ranges.end())
407  {
408  // the two ranges do not overlap or we are at the end of one of the
409  // ranges
410  if (r2 == other.ranges.end() ||
411  (r1 != ranges.end() && r1->end < (r2->begin + offset)))
412  {
413  new_ranges.push_back(*r1);
414  ++r1;
415  }
416  else if (r1 == ranges.end() || (r2->end + offset) < r1->begin)
417  {
418  new_ranges.emplace_back(r2->begin + offset, r2->end + offset);
419  ++r2;
420  }
421  else
422  {
423  // ok, we do overlap, so just take the combination of the current
424  // range (do not bother to merge with subsequent ranges)
425  Range next(std::min(r1->begin, r2->begin + offset),
426  std::max(r1->end, r2->end + offset));
427  new_ranges.push_back(next);
428  ++r1;
429  ++r2;
430  }
431  }
432  ranges.swap(new_ranges);
433 
434  is_compressed = false;
435  compress();
436 }
437 
438 
439 
440 void
441 IndexSet::write(std::ostream &out) const
442 {
443  compress();
444  out << size() << " ";
445  out << ranges.size() << std::endl;
446  std::vector<Range>::const_iterator r = ranges.begin();
447  for (; r != ranges.end(); ++r)
448  {
449  out << r->begin << " " << r->end << std::endl;
450  }
451 }
452 
453 
454 
455 void
456 IndexSet::read(std::istream &in)
457 {
458  AssertThrow(in, ExcIO());
459 
460  size_type s;
461  unsigned int n_ranges;
462 
463  in >> s >> n_ranges;
464  ranges.clear();
465  set_size(s);
466  for (unsigned int i = 0; i < n_ranges; ++i)
467  {
468  AssertThrow(in, ExcIO());
469 
470  size_type b, e;
471  in >> b >> e;
472  add_range(b, e);
473  }
474 }
475 
476 
477 void
478 IndexSet::block_write(std::ostream &out) const
479 {
480  AssertThrow(out, ExcIO());
481  out.write(reinterpret_cast<const char *>(&index_space_size),
482  sizeof(index_space_size));
483  std::size_t n_ranges = ranges.size();
484  out.write(reinterpret_cast<const char *>(&n_ranges), sizeof(n_ranges));
485  if (ranges.empty() == false)
486  out.write(reinterpret_cast<const char *>(&*ranges.begin()),
487  ranges.size() * sizeof(Range));
488  AssertThrow(out, ExcIO());
489 }
490 
491 void
492 IndexSet::block_read(std::istream &in)
493 {
494  size_type size;
495  std::size_t n_ranges;
496  in.read(reinterpret_cast<char *>(&size), sizeof(size));
497  in.read(reinterpret_cast<char *>(&n_ranges), sizeof(n_ranges));
498  // we have to clear ranges first
499  ranges.clear();
500  set_size(size);
501  ranges.resize(n_ranges, Range(0, 0));
502  if (n_ranges != 0u)
503  in.read(reinterpret_cast<char *>(&*ranges.begin()),
504  ranges.size() * sizeof(Range));
505 
506  do_compress(); // needed so that largest_range can be recomputed
507 }
508 
509 
510 
511 void
512 IndexSet::fill_index_vector(std::vector<size_type> &indices) const
513 {
514  compress();
515 
516  indices.clear();
517  indices.reserve(n_elements());
518 
519  for (const auto &range : ranges)
520  for (size_type entry = range.begin; entry < range.end; ++entry)
521  indices.push_back(entry);
522 
523  Assert(indices.size() == n_elements(), ExcInternalError());
524 }
525 
526 
527 
528 #ifdef DEAL_II_WITH_TRILINOS
529 # ifdef DEAL_II_TRILINOS_WITH_TPETRA
530 
531 Tpetra::Map<int, types::global_dof_index>
533  const bool overlapping) const
534 {
535  compress();
536  (void)communicator;
537 
538 # ifdef DEBUG
539  if (!overlapping)
540  {
542  Utilities::MPI::sum(n_elements(), communicator);
543  Assert(n_global_elements == size(),
544  ExcMessage("You are trying to create an Tpetra::Map object "
545  "that partitions elements of an index set "
546  "between processors. However, the union of the "
547  "index sets on different processors does not "
548  "contain all indices exactly once: the sum of "
549  "the number of entries the various processors "
550  "want to store locally is " +
551  std::to_string(n_global_elements) +
552  " whereas the total size of the object to be "
553  "allocated is " +
554  std::to_string(size()) +
555  ". In other words, there are "
556  "either indices that are not spoken for "
557  "by any processor, or there are indices that are "
558  "claimed by multiple processors."));
559  }
560 # endif
561 
562  // Find out if the IndexSet is ascending and 1:1. This corresponds to a
563  // linear Tpetra::Map. Overlapping IndexSets are never 1:1.
564  const bool linear =
565  overlapping ? false : is_ascending_and_one_to_one(communicator);
566  if (linear)
567  return Tpetra::Map<int, types::global_dof_index>(
568  size(),
569  n_elements(),
570  0,
571 # ifdef DEAL_II_WITH_MPI
572  Teuchos::rcp(new Teuchos::MpiComm<int>(communicator))
573 # else
574  Teuchos::rcp(new Teuchos::Comm<int>())
575 # endif
576  );
577  else
578  {
579  std::vector<size_type> indices;
580  fill_index_vector(indices);
581  std::vector<types::global_dof_index> int_indices(indices.size());
582  std::copy(indices.begin(), indices.end(), int_indices.begin());
583  const Teuchos::ArrayView<types::global_dof_index> arr_view(int_indices);
584  return Tpetra::Map<int, types::global_dof_index>(
585  size(),
586  arr_view,
587  0,
588 # ifdef DEAL_II_WITH_MPI
589  Teuchos::rcp(new Teuchos::MpiComm<int>(communicator))
590 # else
591  Teuchos::rcp(new Teuchos::Comm<int>())
592 # endif
593  );
594  }
595 }
596 # endif
597 
598 
599 
600 Epetra_Map
602  const bool overlapping) const
603 {
604  compress();
605  (void)communicator;
606 
607 # ifdef DEBUG
608  if (!overlapping)
609  {
611  Utilities::MPI::sum(n_elements(), communicator);
612  Assert(n_global_elements == size(),
613  ExcMessage("You are trying to create an Epetra_Map object "
614  "that partitions elements of an index set "
615  "between processors. However, the union of the "
616  "index sets on different processors does not "
617  "contain all indices exactly once: the sum of "
618  "the number of entries the various processors "
619  "want to store locally is " +
620  std::to_string(n_global_elements) +
621  " whereas the total size of the object to be "
622  "allocated is " +
623  std::to_string(size()) +
624  ". In other words, there are "
625  "either indices that are not spoken for "
626  "by any processor, or there are indices that are "
627  "claimed by multiple processors."));
628  }
629 # endif
630 
631  // Find out if the IndexSet is ascending and 1:1. This corresponds to a
632  // linear EpetraMap. Overlapping IndexSets are never 1:1.
633  const bool linear =
634  overlapping ? false : is_ascending_and_one_to_one(communicator);
635 
636  if (linear)
637  return Epetra_Map(TrilinosWrappers::types::int_type(size()),
639  0,
640 # ifdef DEAL_II_WITH_MPI
641  Epetra_MpiComm(communicator)
642 # else
643  Epetra_SerialComm()
644 # endif
645  );
646  else
647  {
648  std::vector<size_type> indices;
649  fill_index_vector(indices);
650  return Epetra_Map(
653  (n_elements() > 0 ?
654  reinterpret_cast<TrilinosWrappers::types::int_type *>(
655  indices.data()) :
656  nullptr),
657  0,
658 # ifdef DEAL_II_WITH_MPI
659  Epetra_MpiComm(communicator)
660 # else
661  Epetra_SerialComm()
662 # endif
663  );
664  }
665 }
666 #endif
667 
668 
669 
670 bool
672 {
673  // If the sum of local elements does not add up to the total size,
674  // the IndexSet can't be complete.
676  Utilities::MPI::sum(n_elements(), communicator);
677  if (n_global_elements != size())
678  return false;
679 
680  if (n_global_elements == 0)
681  return true;
682 
683 #ifdef DEAL_II_WITH_MPI
684  // Non-contiguous IndexSets can't be linear.
685  const bool all_contiguous =
686  (Utilities::MPI::min(is_contiguous() ? 1 : 0, communicator) == 1);
687  if (!all_contiguous)
688  return false;
689 
690  bool is_globally_ascending = true;
691  // we know that there is only one interval
692  types::global_dof_index first_local_dof = (n_elements() > 0) ?
693  *(begin_intervals()->begin()) :
695 
696  const unsigned int my_rank = Utilities::MPI::this_mpi_process(communicator);
697  const std::vector<types::global_dof_index> global_dofs =
698  Utilities::MPI::gather(communicator, first_local_dof, 0);
699 
700  if (my_rank == 0)
701  {
702  // find out if the received std::vector is ascending
703  types::global_dof_index index = 0;
704  while (global_dofs[index] == numbers::invalid_dof_index)
705  ++index;
706  types::global_dof_index old_dof = global_dofs[index++];
707  for (; index < global_dofs.size(); ++index)
708  {
709  const types::global_dof_index new_dof = global_dofs[index];
710  if (new_dof != numbers::invalid_dof_index)
711  {
712  if (new_dof <= old_dof)
713  {
714  is_globally_ascending = false;
715  break;
716  }
717  else
718  old_dof = new_dof;
719  }
720  }
721  }
722 
723  // now broadcast the result
724  int is_ascending = is_globally_ascending ? 1 : 0;
725  int ierr = MPI_Bcast(&is_ascending, 1, MPI_INT, 0, communicator);
726  AssertThrowMPI(ierr);
727 
728  return (is_ascending == 1);
729 #else
730  return true;
731 #endif // DEAL_II_WITH_MPI
732 }
733 
734 
735 
736 std::size_t
738 {
742  sizeof(compress_mutex));
743 }
744 
745 
746 
static const unsigned int invalid_unsigned_int
Definition: types.h:196
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1655
ElementIterator end() const
Definition: index_set.h:1585
types::global_dof_index size_type
Definition: cuda_kernels.h:45
static ::ExceptionBase & ExcIO()
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
void block_read(std::istream &in)
Definition: index_set.cc:492
#define AssertIndexRange(index, range)
Definition: exceptions.h:1720
void add_indices(const ForwardIterator &begin, const ForwardIterator &end)
Definition: index_set.h:1708
#define AssertThrow(cond, exc)
Definition: exceptions.h:1571
types::global_dof_index size_type
Definition: index_set.h:83
std::size_t memory_consumption() const
Definition: index_set.cc:737
size_type index_space_size
Definition: index_set.h:972
Tpetra::Map< int, types::global_dof_index > make_tpetra_map(const MPI_Comm &communicator=MPI_COMM_WORLD, const bool overlapping=false) const
Definition: index_set.cc:532
bool is_ascending_and_one_to_one(const MPI_Comm &communicator) const
Definition: index_set.cc:671
Epetra_Map make_trilinos_map(const MPI_Comm &communicator=MPI_COMM_WORLD, const bool overlapping=false) const
Definition: index_set.cc:601
size_type size() const
Definition: index_set.h:1639
static ::ExceptionBase & ExcMessage(std::string arg1)
void block_write(std::ostream &out) const
Definition: index_set.cc:478
std::vector< Range > ranges
Definition: index_set.h:956
std::string compress(const std::string &input)
Definition: utilities.cc:392
T sum(const T &t, const MPI_Comm &mpi_communicator)
void partition(const SparsityPattern &sparsity_pattern, const unsigned int n_partitions, std::vector< unsigned int > &partition_indices, const Partitioner partitioner=Partitioner::metis)
void subtract_set(const IndexSet &other)
Definition: index_set.cc:266
#define Assert(cond, exc)
Definition: exceptions.h:1461
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
void do_compress() const
Definition: index_set.cc:98
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:407
std::vector< T > gather(const MPI_Comm &comm, const T &object_to_send, const unsigned int root_process=0)
std::string to_string(const T &t)
Definition: patterns.h:2330
std::enable_if< IsBlockVector< VectorType >::value, unsigned int >::type n_blocks(const VectorType &vector)
Definition: operators.h:50
IndexSet tensor_product(const IndexSet &other) const
Definition: index_set.cc:335
size_type pop_back()
Definition: index_set.cc:347
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
void fill_index_vector(std::vector< size_type > &indices) const
Definition: index_set.cc:512
void add_range(const size_type begin, const size_type end)
Definition: index_set.h:1678
TrilinosWrappers::types::int_type n_global_elements(const Epetra_BlockMap &map)
IndexSet operator &(const IndexSet &is) const
void set_size(const size_type size)
Definition: index_set.h:1627
std::vector< IndexSet > split_by_block(const std::vector< types::global_dof_index > &n_indices_per_block) const
Definition: index_set.cc:239
IndexSet get_view(const size_type begin, const size_type end) const
Definition: index_set.cc:212
void compress() const
Definition: index_set.h:1647
IntervalIterator begin_intervals() const
Definition: index_set.h:1594
#define AssertThrowMPI(error_code)
Definition: exceptions.h:1778
bool is_contiguous() const
Definition: index_set.h:1820
Threads::Mutex compress_mutex
Definition: index_set.h:989
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:406
T min(const T &t, const MPI_Comm &mpi_communicator)
size_type pop_front()
Definition: index_set.cc:365
bool is_compressed
Definition: index_set.h:966
ElementIterator begin() const
Definition: index_set.h:1078
unsigned int this_mpi_process(const MPI_Comm &mpi_communicator)
Definition: mpi.cc:128
void write(std::ostream &out) const
Definition: index_set.cc:441
const types::global_dof_index invalid_dof_index
Definition: types.h:211
ElementIterator begin() const
Definition: index_set.h:1523
void read(std::istream &in)
Definition: index_set.cc:456
size_type n_elements() const
Definition: index_set.h:1837
void copy(const T *begin, const T *end, U *dest)
bool is_empty() const
Definition: index_set.h:1829
size_type largest_range
Definition: index_set.h:983
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
size_type nth_index_in_set
Definition: index_set.h:884
static ::ExceptionBase & ExcInternalError()