16 #include <deal.II/base/memory_consumption.h> 17 #include <deal.II/base/mpi.h> 18 #include <deal.II/base/index_set.h> 22 #ifdef DEAL_II_WITH_TRILINOS 23 # ifdef DEAL_II_WITH_MPI 24 # include <Epetra_MpiComm.h> 26 # include <Epetra_SerialComm.h> 27 # include <Epetra_Map.h> 30 DEAL_II_NAMESPACE_OPEN
34 #ifdef DEAL_II_WITH_TRILINOS 39 #ifdef DEAL_II_WITH_64BIT_INDICES 44 index_space_size (map.NumGlobalElements64()),
45 largest_range (
numbers::invalid_unsigned_int)
52 const size_type n_indices = map.NumMyElements();
54 add_indices(indices, indices+n_indices);
66 index_space_size (map.NumGlobalElements()),
67 largest_range (
numbers::invalid_unsigned_int)
74 const size_type n_indices = map.NumMyElements();
75 unsigned int *indices = (
unsigned int *)map.MyGlobalElements();
83 #endif // ifdef DEAL_II_WITH_TRILINOS 98 ExcIndexRangeType<size_type> (
begin, 0,
end));
107 ranges.push_back(new_range);
133 std::vector<Range>::iterator store =
ranges.begin();
134 for (std::vector<Range>::iterator i =
ranges.begin();
137 std::vector<Range>::iterator
145 while (next !=
ranges.end() &&
146 (next->begin <= last_index))
148 last_index = std::max (last_index, next->end);
154 *store =
Range(first_index, last_index);
158 if (store !=
ranges.end())
160 std::vector<Range> new_ranges(
ranges.begin(), store);
165 size_type next_index = 0, largest_range_size = 0;
166 for (std::vector<Range>::iterator i =
ranges.begin(); i !=
ranges.end();
171 i->nth_index_in_set = next_index;
172 next_index += (i->end - i->begin);
173 if (i->end - i->begin > largest_range_size)
175 largest_range_size = i->end - i->begin;
197 std::vector<Range>::const_iterator r1 =
ranges.begin(),
201 while ((r1 !=
ranges.end())
207 if (r1->end <= r2->begin)
209 else if (r2->end <= r1->begin)
214 Assert (((r1->begin <= r2->begin) &&
215 (r1->end > r2->begin))
217 ((r2->begin <= r1->begin) &&
218 (r2->end > r1->begin)),
222 result.add_range (std::max (r1->begin,
230 if (r1->end <= r2->end)
248 ExcMessage (
"End index needs to be larger or equal to begin index!"));
250 ExcMessage (
"Given range exceeds index set dimension"));
253 std::vector<Range>::const_iterator r1 =
ranges.begin();
255 while (r1 !=
ranges.end())
257 if ((r1->end >
begin)
265 else if (r1->begin >=
end)
287 std::vector<Range> new_ranges;
289 std::vector<Range>::iterator own_it =
ranges.begin();
290 std::vector<Range>::iterator other_it = other.
ranges.begin();
292 while (own_it !=
ranges.end() && other_it != other.
ranges.end())
295 if (own_it->end <= other_it->begin)
301 if (own_it->begin >= other_it->end)
309 if (own_it->begin < other_it->begin)
311 Range r(own_it->begin, other_it->begin);
312 r.nth_index_in_set = 0;
313 new_ranges.push_back(r);
318 own_it->begin = other_it->end;
319 if (own_it->begin > own_it->end)
321 own_it->begin = own_it->end;
331 for (std::vector<Range>::iterator it =
ranges.begin();
334 if (it->begin >= it->end)
341 const std::vector<Range>::iterator
end = new_ranges.end();
342 for (std::vector<Range>::iterator it = new_ranges.begin();
355 ExcMessage(
"pop_back() failed, because this IndexSet contains no entries."));
372 ExcMessage(
"pop_front() failed, because this IndexSet contains no entries."));
391 const unsigned int offset)
393 if ((
this == &other) && (offset == 0))
398 ExcIndexRangeType<size_type> (other.
ranges.back().end-1,
404 std::vector<Range>::const_iterator r1 =
ranges.begin(),
405 r2 = other.
ranges.begin();
407 std::vector<Range> new_ranges;
414 if (r2 == other.
ranges.end() ||
415 (r1 !=
ranges.end() && r1->end < (r2->begin+offset)))
417 new_ranges.push_back(*r1);
420 else if (r1 ==
ranges.end() || (r2->end+offset) < r1->begin)
422 new_ranges.push_back(
Range(r2->begin+offset,r2->end+offset));
429 Range next(std::min(r1->begin, r2->begin+offset),
430 std::max(r1->end, r2->end+offset));
431 new_ranges.push_back(next);
448 out <<
size() <<
" ";
449 out <<
ranges.size() << std::endl;
450 std::vector<Range>::const_iterator r =
ranges.begin();
451 for ( ; r!=
ranges.end(); ++r)
453 out << r->begin <<
" " << r->end << std::endl;
463 unsigned int numranges;
465 in >> s >> numranges;
468 for (
unsigned int i=0; i<numranges; ++i)
483 size_t n_ranges =
ranges.size();
484 out.write(reinterpret_cast<const char *>(&n_ranges),
486 if (
ranges.empty() ==
false)
487 out.write (reinterpret_cast<const char *>(&*
ranges.begin()),
497 in.read(reinterpret_cast<char *>(&
size),
sizeof(
size));
498 in.read(reinterpret_cast<char *>(&n_ranges),
sizeof(n_ranges));
504 in.read(reinterpret_cast<char *>(&*
ranges.begin()),
519 for (std::vector<Range>::iterator it =
ranges.begin();
523 indices.push_back (i);
532 #ifdef DEAL_II_WITH_TRILINOS 536 const bool overlapping)
const 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 " 555 " whereas the total size of the object to be " 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."));
570 return Epetra_Map (TrilinosWrappers::types::int_type(
size()),
571 TrilinosWrappers::types::int_type(
n_elements()),
573 #ifdef DEAL_II_WITH_MPI
574 Epetra_MpiComm(communicator)
581 std::vector<size_type> indices;
583 return Epetra_Map (TrilinosWrappers::types::int_type(-1),
584 TrilinosWrappers::types::int_type(
n_elements()),
587 reinterpret_cast<TrilinosWrappers::types::int_type *>(&indices[0])
591 #ifdef DEAL_II_WITH_MPI
592 Epetra_MpiComm(communicator)
610 if (n_global_elements !=
size())
613 #ifdef DEAL_II_WITH_MPI 619 bool is_globally_ascending =
true;
628 const unsigned int gather_size = (my_rank==0)?n_ranks:1;
629 std::vector<types::global_dof_index> global_dofs(gather_size);
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,
643 old_dof = global_dofs[index++];
644 for (; index<global_dofs.size(); ++index)
649 if (new_dof <= old_dof)
651 is_globally_ascending =
false;
661 int is_ascending = is_globally_ascending ? 1 : 0;
662 ierr = MPI_Bcast(&is_ascending, 1, MPI_INT, 0, communicator);
665 return (is_ascending==1);
668 #endif //DEAL_II_WITH_MPI 684 DEAL_II_NAMESPACE_CLOSE
Iterator lower_bound(Iterator first, Iterator last, const T &val)
types::global_dof_index size_type
ElementIterator end() const
IndexSet operator&(const IndexSet &is) const
static ::ExceptionBase & ExcIO()
void block_read(std::istream &in)
void add_indices(const ForwardIterator &begin, const ForwardIterator &end)
#define AssertThrow(cond, exc)
std::string to_string(const number value, const unsigned int digits=numbers::invalid_unsigned_int)
std::size_t memory_consumption() const
size_type index_space_size
bool is_ascending_and_one_to_one(const MPI_Comm &communicator) const
Epetra_Map make_trilinos_map(const MPI_Comm &communicator=MPI_COMM_WORLD, const bool overlapping=false) const
static ::ExceptionBase & ExcMessage(std::string arg1)
void block_write(std::ostream &out) const
std::vector< Range > ranges
T sum(const T &t, const MPI_Comm &mpi_communicator)
unsigned int global_dof_index
void subtract_set(const IndexSet &other)
#define Assert(cond, exc)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
unsigned int n_mpi_processes(const MPI_Comm &mpi_communicator)
void fill_index_vector(std::vector< size_type > &indices) const
void add_range(const size_type begin, const size_type end)
void set_size(const size_type size)
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
IntervalIterator begin_intervals() const
#define AssertThrowMPI(error_code)
bool is_contiguous() const
Threads::Mutex compress_mutex
T min(const T &t, const MPI_Comm &mpi_communicator)
ElementIterator begin() const
unsigned int this_mpi_process(const MPI_Comm &mpi_communicator)
void write(std::ostream &out) const
const types::global_dof_index invalid_dof_index
types::global_dof_index size_type
ElementIterator begin() const
void read(std::istream &in)
size_type n_elements() const
static ::ExceptionBase & ExcInternalError()