Reference documentation for deal.II version 8.5.1
index_set.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2005 - 2017 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 at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #include <deal.II/base/memory_consumption.h>
17 #include <deal.II/base/mpi.h>
18 #include <deal.II/base/index_set.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_SerialComm.h>
27 # include <Epetra_Map.h>
28 #endif
29 
30 DEAL_II_NAMESPACE_OPEN
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_Map &map)
42  :
43  is_compressed (true),
44  index_space_size (map.NumGlobalElements64()),
45  largest_range (numbers::invalid_unsigned_int)
46 {
47  // For a contiguous map, we do not need to go through the whole data...
48  if (map.LinearMap())
49  add_range(size_type(map.MinMyGID64()), size_type(map.MaxMyGID64()+1));
50  else
51  {
52  const size_type n_indices = map.NumMyElements();
53  size_type *indices = (size_type *)map.MyGlobalElements64();
54  add_indices(indices, indices+n_indices);
55  }
56  compress();
57 }
58 
59 #else
60 
61 // this is the standard 32-bit implementation
62 
63 IndexSet::IndexSet (const Epetra_Map &map)
64  :
65  is_compressed (true),
66  index_space_size (map.NumGlobalElements()),
67  largest_range (numbers::invalid_unsigned_int)
68 {
69  // For a contiguous map, we do not need to go through the whole data...
70  if (map.LinearMap())
71  add_range(size_type(map.MinMyGID()), size_type(map.MaxMyGID()+1));
72  else
73  {
74  const size_type n_indices = map.NumMyElements();
75  unsigned int *indices = (unsigned int *)map.MyGlobalElements();
76  add_indices(indices, indices+n_indices);
77  }
78  compress();
79 }
80 
81 #endif
82 
83 #endif // ifdef DEAL_II_WITH_TRILINOS
84 
85 
86 
87 void
89  const size_type end)
90 {
92  ||
94  ExcIndexRangeType<size_type> (begin, 0, index_space_size));
96  ExcIndexRangeType<size_type> (end, 0, index_space_size+1));
97  Assert (begin <= end,
98  ExcIndexRangeType<size_type> (begin, 0, end));
99 
100  if (begin != end)
101  {
102  const Range new_range(begin,end);
103 
104  // the new index might be larger than the last index present in the
105  // ranges. Then we can skip the binary search
106  if (ranges.size() == 0 || begin > ranges.back().end)
107  ranges.push_back(new_range);
108  else
109  ranges.insert (Utilities::lower_bound (ranges.begin(),
110  ranges.end(),
111  new_range),
112  new_range);
113  is_compressed = false;
114  }
115 }
116 
117 
118 
119 void
121 {
122  // we will, in the following, modify mutable variables. this can only
123  // work in multithreaded applications if we lock the data structures
124  // via a mutex, so that users can call 'const' functions from threads
125  // in parallel (and these 'const' functions can then call compress()
126  // which itself calls the current function)
128 
129  // see if any of the contiguous ranges can be merged. do not use
130  // std::vector::erase in-place as it is quadratic in the number of
131  // ranges. since the ranges are sorted by their first index, determining
132  // overlap isn't all that hard
133  std::vector<Range>::iterator store = ranges.begin();
134  for (std::vector<Range>::iterator i = ranges.begin();
135  i != ranges.end(); )
136  {
137  std::vector<Range>::iterator
138  next = i;
139  ++next;
140 
141  size_type first_index = i->begin;
142  size_type last_index = i->end;
143 
144  // see if we can merge any of the following ranges
145  while (next != ranges.end() &&
146  (next->begin <= last_index))
147  {
148  last_index = std::max (last_index, next->end);
149  ++next;
150  }
151  i = next;
152 
153  // store the new range in the slot we last occupied
154  *store = Range(first_index, last_index);
155  ++store;
156  }
157  // use a compact array with exactly the right amount of storage
158  if (store != ranges.end())
159  {
160  std::vector<Range> new_ranges(ranges.begin(), store);
161  ranges.swap(new_ranges);
162  }
163 
164  // now compute indices within set and the range with most elements
165  size_type next_index = 0, largest_range_size = 0;
166  for (std::vector<Range>::iterator i = ranges.begin(); i != ranges.end();
167  ++i)
168  {
169  Assert(i->begin < i->end, ExcInternalError());
170 
171  i->nth_index_in_set = next_index;
172  next_index += (i->end - i->begin);
173  if (i->end - i->begin > largest_range_size)
174  {
175  largest_range_size = i->end - i->begin;
176  largest_range = i - ranges.begin();
177  }
178  }
179  is_compressed = true;
180 
181  // check that next_index is correct. needs to be after the previous
182  // statement because we otherwise will get into an endless loop
183  Assert (next_index == n_elements(), ExcInternalError());
184 }
185 
186 
187 
188 IndexSet
189 IndexSet::operator & (const IndexSet &is) const
190 {
191  Assert (size() == is.size(),
192  ExcDimensionMismatch (size(), is.size()));
193 
194  compress ();
195  is.compress ();
196 
197  std::vector<Range>::const_iterator r1 = ranges.begin(),
198  r2 = is.ranges.begin();
199  IndexSet result (size());
200 
201  while ((r1 != ranges.end())
202  &&
203  (r2 != is.ranges.end()))
204  {
205  // if r1 and r2 do not overlap at all, then move the pointer that sits
206  // to the left of the other up by one
207  if (r1->end <= r2->begin)
208  ++r1;
209  else if (r2->end <= r1->begin)
210  ++r2;
211  else
212  {
213  // the ranges must overlap somehow
214  Assert (((r1->begin <= r2->begin) &&
215  (r1->end > r2->begin))
216  ||
217  ((r2->begin <= r1->begin) &&
218  (r2->end > r1->begin)),
219  ExcInternalError());
220 
221  // add the overlapping range to the result
222  result.add_range (std::max (r1->begin,
223  r2->begin),
224  std::min (r1->end,
225  r2->end));
226 
227  // now move that iterator that ends earlier one up. note that it has
228  // to be this one because a subsequent range may still have a chance
229  // of overlapping with the range that ends later
230  if (r1->end <= r2->end)
231  ++r1;
232  else
233  ++r2;
234  }
235  }
236 
237  result.compress ();
238  return result;
239 }
240 
241 
242 
243 IndexSet
245  const size_type end) const
246 {
247  Assert (begin <= end,
248  ExcMessage ("End index needs to be larger or equal to begin index!"));
249  Assert (end <= size(),
250  ExcMessage ("Given range exceeds index set dimension"));
251 
252  IndexSet result (end-begin);
253  std::vector<Range>::const_iterator r1 = ranges.begin();
254 
255  while (r1 != ranges.end())
256  {
257  if ((r1->end > begin)
258  &&
259  (r1->begin < end))
260  {
261  result.add_range (std::max(r1->begin, begin)-begin,
262  std::min(r1->end, end)-begin);
263 
264  }
265  else if (r1->begin >= end)
266  break;
267 
268  ++r1;
269  }
270 
271  result.compress();
272  return result;
273 }
274 
275 
276 
277 void
279 {
280  compress();
281  other.compress();
282  is_compressed = false;
283 
284 
285  // we save new ranges to be added to our IndexSet in an temporary vector and
286  // add all of them in one go at the end.
287  std::vector<Range> new_ranges;
288 
289  std::vector<Range>::iterator own_it = ranges.begin();
290  std::vector<Range>::iterator other_it = other.ranges.begin();
291 
292  while (own_it != ranges.end() && other_it != other.ranges.end())
293  {
294  //advance own iterator until we get an overlap
295  if (own_it->end <= other_it->begin)
296  {
297  ++own_it;
298  continue;
299  }
300  //we are done with other_it, so advance
301  if (own_it->begin >= other_it->end)
302  {
303  ++other_it;
304  continue;
305  }
306 
307  //Now own_it and other_it overlap. First save the part of own_it that
308  //is before other_it (if not empty).
309  if (own_it->begin < other_it->begin)
310  {
311  Range r(own_it->begin, other_it->begin);
312  r.nth_index_in_set = 0; //fix warning of unused variable
313  new_ranges.push_back(r);
314  }
315  // change own_it to the sub range behind other_it. Do not delete own_it
316  // in any case. As removal would invalidate iterators, we just shrink
317  // the range to an empty one.
318  own_it->begin = other_it->end;
319  if (own_it->begin > own_it->end)
320  {
321  own_it->begin = own_it->end;
322  ++own_it;
323  }
324 
325  // continue without advancing iterators, the right one will be advanced
326  // next.
327  }
328 
329  // Now delete all empty ranges we might
330  // have created.
331  for (std::vector<Range>::iterator it = ranges.begin();
332  it != ranges.end(); )
333  {
334  if (it->begin >= it->end)
335  it = ranges.erase(it);
336  else
337  ++it;
338  }
339 
340  // done, now add the temporary ranges
341  const std::vector<Range>::iterator end = new_ranges.end();
342  for (std::vector<Range>::iterator it = new_ranges.begin();
343  it != end; ++it)
344  add_range(it->begin, it->end);
345 
346  compress();
347 }
348 
349 
350 
353 {
354  Assert(is_empty() == false,
355  ExcMessage("pop_back() failed, because this IndexSet contains no entries."));
356 
357  const size_type index = ranges.back().end-1;
358  --ranges.back().end;
359 
360  if (ranges.back().begin == ranges.back().end)
361  ranges.pop_back();
362 
363  return index;
364 }
365 
366 
367 
370 {
371  Assert(is_empty() == false,
372  ExcMessage("pop_front() failed, because this IndexSet contains no entries."));
373 
374  const size_type index = ranges.front().begin;
375  ++ranges.front().begin;
376 
377  if (ranges.front().begin == ranges.front().end)
378  ranges.erase(ranges.begin());
379 
380  // We have to set this in any case, because nth_index_in_set is no longer
381  // up to date for all but the first range
382  is_compressed = false;
383 
384  return index;
385 }
386 
387 
388 
389 void
391  const unsigned int offset)
392 {
393  if ((this == &other) && (offset == 0))
394  return;
395 
396  Assert (other.ranges.size() == 0
397  || other.ranges.back().end-1 < index_space_size,
398  ExcIndexRangeType<size_type> (other.ranges.back().end-1,
399  0, index_space_size));
400 
401  compress();
402  other.compress();
403 
404  std::vector<Range>::const_iterator r1 = ranges.begin(),
405  r2 = other.ranges.begin();
406 
407  std::vector<Range> new_ranges;
408  // just get the start and end of the ranges right in this method, everything
409  // else will be done in compress()
410  while (r1 != ranges.end() || r2 != other.ranges.end())
411  {
412  // the two ranges do not overlap or we are at the end of one of the
413  // ranges
414  if (r2 == other.ranges.end() ||
415  (r1 != ranges.end() && r1->end < (r2->begin+offset)))
416  {
417  new_ranges.push_back(*r1);
418  ++r1;
419  }
420  else if (r1 == ranges.end() || (r2->end+offset) < r1->begin)
421  {
422  new_ranges.push_back(Range(r2->begin+offset,r2->end+offset));
423  ++r2;
424  }
425  else
426  {
427  // ok, we do overlap, so just take the combination of the current
428  // range (do not bother to merge with subsequent ranges)
429  Range next(std::min(r1->begin, r2->begin+offset),
430  std::max(r1->end, r2->end+offset));
431  new_ranges.push_back(next);
432  ++r1;
433  ++r2;
434  }
435  }
436  ranges.swap(new_ranges);
437 
438  is_compressed = false;
439  compress();
440 }
441 
442 
443 
444 void
445 IndexSet::write(std::ostream &out) const
446 {
447  compress();
448  out << size() << " ";
449  out << ranges.size() << std::endl;
450  std::vector<Range>::const_iterator r = ranges.begin();
451  for ( ; r!=ranges.end(); ++r)
452  {
453  out << r->begin << " " << r->end << std::endl;
454  }
455 }
456 
457 
458 
459 void
460 IndexSet::read(std::istream &in)
461 {
462  size_type s;
463  unsigned int numranges;
464 
465  in >> s >> numranges;
466  ranges.clear();
467  set_size(s);
468  for (unsigned int i=0; i<numranges; ++i)
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  size_t n_ranges = ranges.size();
484  out.write(reinterpret_cast<const char *>(&n_ranges),
485  sizeof(n_ranges));
486  if (ranges.empty() == false)
487  out.write (reinterpret_cast<const char *>(&*ranges.begin()),
488  ranges.size() * sizeof(Range));
489  AssertThrow (out, ExcIO());
490 }
491 
492 void
493 IndexSet::block_read(std::istream &in)
494 {
495  size_type size;
496  size_t n_ranges;
497  in.read(reinterpret_cast<char *>(&size), sizeof(size));
498  in.read(reinterpret_cast<char *>(&n_ranges), sizeof(n_ranges));
499  // we have to clear ranges first
500  ranges.clear();
501  set_size(size);
502  ranges.resize(n_ranges, Range(0,0));
503  if (n_ranges)
504  in.read(reinterpret_cast<char *>(&*ranges.begin()),
505  ranges.size() * sizeof(Range));
506 
507  do_compress(); // needed so that largest_range can be recomputed
508 }
509 
510 
511 
512 void 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 (std::vector<Range>::iterator it = ranges.begin();
520  it != ranges.end();
521  ++it)
522  for (size_type i=it->begin; i<it->end; ++i)
523  indices.push_back (i);
524 
525  Assert (indices.size() == n_elements(), ExcInternalError());
526 }
527 
528 
529 
530 
531 
532 #ifdef DEAL_II_WITH_TRILINOS
533 
534 Epetra_Map
535 IndexSet::make_trilinos_map (const MPI_Comm &communicator,
536  const bool overlapping) const
537 {
538  compress ();
539  (void)communicator;
540 
541 #ifdef DEBUG
542  if (!overlapping)
543  {
544  const size_type n_global_elements
545  = Utilities::MPI::sum (n_elements(), communicator);
546  Assert (n_global_elements == size(),
547  ExcMessage ("You are trying to create an Epetra_Map object "
548  "that partitions elements of an index set "
549  "between processors. However, the union of the "
550  "index sets on different processors does not "
551  "contain all indices exactly once: the sum of "
552  "the number of entries the various processors "
553  "want to store locally is "
554  + Utilities::to_string (n_global_elements) +
555  " whereas the total size of the object to be "
556  "allocated is "
557  + Utilities::to_string (size()) +
558  ". In other words, there are "
559  "either indices that are not spoken for "
560  "by any processor, or there are indices that are "
561  "claimed by multiple processors."));
562  }
563 #endif
564 
565  // Find out if the IndexSet is ascending and 1:1. This correspinds to a
566  // linear EpetraMap. Overlapping IndexSets are never 1:1.
567  const bool linear = overlapping ? false : is_ascending_and_one_to_one(communicator);
568 
569  if (linear)
570  return Epetra_Map (TrilinosWrappers::types::int_type(size()),
571  TrilinosWrappers::types::int_type(n_elements()),
572  0,
573 #ifdef DEAL_II_WITH_MPI
574  Epetra_MpiComm(communicator)
575 #else
576  Epetra_SerialComm()
577 #endif
578  );
579  else
580  {
581  std::vector<size_type> indices;
582  fill_index_vector(indices);
583  return Epetra_Map (TrilinosWrappers::types::int_type(-1),
584  TrilinosWrappers::types::int_type(n_elements()),
585  (n_elements() > 0
586  ?
587  reinterpret_cast<TrilinosWrappers::types::int_type *>(&indices[0])
588  :
589  0),
590  0,
591 #ifdef DEAL_II_WITH_MPI
592  Epetra_MpiComm(communicator)
593 #else
594  Epetra_SerialComm()
595 #endif
596  );
597  }
598 }
599 #endif
600 
601 
602 
603 bool
604 IndexSet::is_ascending_and_one_to_one (const MPI_Comm &communicator) const
605 {
606  // If the sum of local elements does not add up to the total size,
607  // the IndexSet can't be complete.
608  const size_type n_global_elements
609  = Utilities::MPI::sum (n_elements(), communicator);
610  if (n_global_elements != size())
611  return false;
612 
613 #ifdef DEAL_II_WITH_MPI
614  // Non-contiguous IndexSets can't be linear.
615  const bool all_contiguous = (Utilities::MPI::min (is_contiguous() ? 1 : 0, communicator) == 1);
616  if (!all_contiguous)
617  return false;
618 
619  bool is_globally_ascending = true;
620  // we know that there is only one interval
621  types::global_dof_index first_local_dof
622  = (n_elements()>0) ? *(begin_intervals()->begin())
624 
625  const unsigned int my_rank = Utilities::MPI::this_mpi_process(communicator);
626  const unsigned int n_ranks = Utilities::MPI::n_mpi_processes(communicator);
627  // first gather all information on process 0
628  const unsigned int gather_size = (my_rank==0)?n_ranks:1;
629  std::vector<types::global_dof_index> global_dofs(gather_size);
630 
631  int ierr = MPI_Gather(&first_local_dof, 1, DEAL_II_DOF_INDEX_MPI_TYPE,
632  &(global_dofs[0]), 1, DEAL_II_DOF_INDEX_MPI_TYPE, 0,
633  communicator);
634  AssertThrowMPI(ierr);
635 
636  if (my_rank == 0)
637  {
638  // find out if the received std::vector is ascending
639  types::global_dof_index old_dof = global_dofs[0];
640  types::global_dof_index index = 0;
641  while (global_dofs[index] == numbers::invalid_dof_index)
642  ++index;
643  old_dof = global_dofs[index++];
644  for (; index<global_dofs.size(); ++index)
645  {
646  const types::global_dof_index new_dof = global_dofs[index];
647  if (new_dof != numbers::invalid_dof_index)
648  {
649  if (new_dof <= old_dof)
650  {
651  is_globally_ascending = false;
652  break;
653  }
654  else
655  old_dof = new_dof;
656  }
657  }
658  }
659 
660  // now broadcast the result
661  int is_ascending = is_globally_ascending ? 1 : 0;
662  ierr = MPI_Bcast(&is_ascending, 1, MPI_INT, 0, communicator);
663  AssertThrowMPI(ierr);
664 
665  return (is_ascending==1);
666 #else
667  return true;
668 #endif //DEAL_II_WITH_MPI
669 }
670 
671 
672 
673 std::size_t
675 {
679  sizeof (compress_mutex));
680 }
681 
682 
683 
684 DEAL_II_NAMESPACE_CLOSE
Iterator lower_bound(Iterator first, Iterator last, const T &val)
Definition: utilities.h:618
ElementIterator end() const
Definition: index_set.h:1305
IndexSet operator&(const IndexSet &is) const
static ::ExceptionBase & ExcIO()
void block_read(std::istream &in)
Definition: index_set.cc:493
void add_indices(const ForwardIterator &begin, const ForwardIterator &end)
Definition: index_set.h:1463
#define AssertThrow(cond, exc)
Definition: exceptions.h:369
std::string to_string(const number value, const unsigned int digits=numbers::invalid_unsigned_int)
Definition: utilities.cc:92
std::size_t memory_consumption() const
Definition: index_set.cc:674
size_type index_space_size
Definition: index_set.h:844
bool is_ascending_and_one_to_one(const MPI_Comm &communicator) const
Definition: index_set.cc:604
Epetra_Map make_trilinos_map(const MPI_Comm &communicator=MPI_COMM_WORLD, const bool overlapping=false) const
Definition: index_set.cc:535
size_type size() const
Definition: index_set.h:1419
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:828
T sum(const T &t, const MPI_Comm &mpi_communicator)
unsigned int global_dof_index
Definition: types.h:88
void subtract_set(const IndexSet &other)
Definition: index_set.cc:278
#define Assert(cond, exc)
Definition: exceptions.h:313
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
void do_compress() const
Definition: index_set.cc:120
unsigned int n_mpi_processes(const MPI_Comm &mpi_communicator)
Definition: mpi.cc:63
size_type pop_back()
Definition: index_set.cc:352
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.cc:88
void set_size(const size_type size)
Definition: index_set.h:1406
std_cxx11::enable_if< std_cxx11::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
IndexSet get_view(const size_type begin, const size_type end) const
Definition: index_set.cc:244
void compress() const
Definition: index_set.h:1428
IntervalIterator begin_intervals() const
Definition: index_set.h:1313
#define AssertThrowMPI(error_code)
Definition: exceptions.h:1211
bool is_contiguous() const
Definition: index_set.h:1542
Threads::Mutex compress_mutex
Definition: index_set.h:861
T min(const T &t, const MPI_Comm &mpi_communicator)
size_type pop_front()
Definition: index_set.cc:369
bool is_compressed
Definition: index_set.h:838
ElementIterator begin() const
Definition: index_set.h:926
unsigned int this_mpi_process(const MPI_Comm &mpi_communicator)
Definition: mpi.cc:73
void write(std::ostream &out) const
Definition: index_set.cc:445
const types::global_dof_index invalid_dof_index
Definition: types.h:184
types::global_dof_index size_type
Definition: index_set.h:73
ElementIterator begin() const
Definition: index_set.h:1246
void read(std::istream &in)
Definition: index_set.cc:460
size_type n_elements() const
Definition: index_set.h:1560
bool is_empty() const
Definition: index_set.h:1551
size_type largest_range
Definition: index_set.h:855
static ::ExceptionBase & ExcInternalError()