16 #include <deal.II/base/memory_consumption.h> 17 #include <deal.II/base/index_set.h> 20 #ifdef DEAL_II_WITH_TRILINOS 21 # ifdef DEAL_II_WITH_MPI 22 # include <Epetra_MpiComm.h> 24 # include <Epetra_SerialComm.h> 25 # include <Epetra_Map.h> 28 DEAL_II_NAMESPACE_OPEN
32 #ifdef DEAL_II_WITH_TRILINOS 37 #ifdef DEAL_II_WITH_64BIT_INDICES 42 index_space_size (map.NumGlobalElements64()),
43 largest_range (
numbers::invalid_unsigned_int)
50 const size_type n_indices = map.NumMyElements();
52 add_indices(indices, indices+n_indices);
64 index_space_size (map.NumGlobalElements()),
65 largest_range (
numbers::invalid_unsigned_int)
72 const size_type n_indices = map.NumMyElements();
73 unsigned int *indices = (
unsigned int *)map.MyGlobalElements();
81 #endif // ifdef DEAL_II_WITH_TRILINOS 96 ExcIndexRangeType<size_type> (begin, 0, end));
100 const Range new_range(begin,end);
105 ranges.push_back(new_range);
124 std::vector<Range>::iterator store =
ranges.begin();
125 for (std::vector<Range>::iterator i =
ranges.begin();
128 std::vector<Range>::iterator
136 while (next !=
ranges.end() &&
137 (next->begin <= last_index))
139 last_index = std::max (last_index, next->end);
145 *store =
Range(first_index, last_index);
149 if (store !=
ranges.end())
151 std::vector<Range> new_ranges(
ranges.begin(), store);
156 size_type next_index = 0, largest_range_size = 0;
157 for (std::vector<Range>::iterator i =
ranges.begin(); i !=
ranges.end();
160 Assert(i->begin < i->end, ExcInternalError());
162 i->nth_index_in_set = next_index;
163 next_index += (i->end - i->begin);
164 if (i->end - i->begin > largest_range_size)
166 largest_range_size = i->end - i->begin;
183 ExcDimensionMismatch (
size(), is.
size()));
188 std::vector<Range>::const_iterator r1 =
ranges.begin(),
192 while ((r1 !=
ranges.end())
198 if (r1->end <= r2->begin)
200 else if (r2->end <= r1->begin)
205 Assert (((r1->begin <= r2->begin) &&
206 (r1->end > r2->begin))
208 ((r2->begin <= r1->begin) &&
209 (r2->end > r1->begin)),
221 if (r1->end <= r2->end)
239 ExcMessage (
"End index needs to be larger or equal to begin index!"));
241 ExcMessage (
"Given range exceeds index set dimension"));
244 std::vector<Range>::const_iterator r1 =
ranges.begin();
246 while (r1 !=
ranges.end())
248 if ((r1->end > begin)
252 result.
add_range (std::max(r1->begin, begin)-begin,
253 std::min(r1->end, end)-begin);
256 else if (r1->begin >= end)
279 std::list<Range> temp_list;
281 std::vector<Range>::iterator own_it =
ranges.begin();
282 std::vector<Range>::iterator other_it = other.
ranges.begin();
284 while (own_it !=
ranges.end() && other_it != other.
ranges.end())
287 if (own_it->end <= other_it->begin)
293 if (own_it->begin >= other_it->end)
301 if (own_it->begin < other_it->begin)
303 Range r(own_it->begin, other_it->begin);
304 r.nth_index_in_set = 0;
305 temp_list.push_back(r);
310 own_it->begin = other_it->end;
311 if (own_it->begin > own_it->end)
313 own_it->begin = own_it->end;
323 for (std::vector<Range>::iterator it =
ranges.begin();
326 if (it->begin >= it->end)
333 for (std::list<Range>::iterator it = temp_list.begin();
334 it != temp_list.end();
345 const unsigned int offset)
347 if ((
this == &other) && (offset == 0))
353 std::vector<Range>::const_iterator r1 =
ranges.begin(),
354 r2 = other.
ranges.begin();
356 std::vector<Range> new_ranges;
363 if (r2 == other.
ranges.end() ||
364 (r1 !=
ranges.end() && r1->end < (r2->begin+offset)))
366 new_ranges.push_back(*r1);
369 else if (r1 ==
ranges.end() || (r2->end+offset) < r1->begin)
371 new_ranges.push_back(
Range(r2->begin+offset,r2->end+offset));
378 Range next(std::min(r1->begin, r2->begin+offset),
379 std::max(r1->end, r2->end+offset));
380 new_ranges.push_back(next);
397 out <<
size() <<
" ";
398 out <<
ranges.size() << std::endl;
399 std::vector<Range>::const_iterator r =
ranges.begin();
400 for ( ; r!=
ranges.end(); ++r)
402 out << r->begin <<
" " << r->end << std::endl;
412 unsigned int numranges;
414 in >> s >> numranges;
417 for (
unsigned int i=0; i<numranges; ++i)
432 size_t n_ranges =
ranges.size();
433 out.write(reinterpret_cast<const char *>(&n_ranges),
435 if (
ranges.empty() ==
false)
436 out.write (reinterpret_cast<const char *>(&*
ranges.begin()),
446 in.read(reinterpret_cast<char *>(&size),
sizeof(size));
447 in.read(reinterpret_cast<char *>(&n_ranges),
sizeof(n_ranges));
453 in.read(reinterpret_cast<char *>(&*
ranges.begin()),
468 for (std::vector<Range>::iterator it =
ranges.begin();
472 indices.push_back (i);
481 #ifdef DEAL_II_WITH_TRILINOS 485 const bool overlapping)
const 495 ExcMessage (
"You are trying to create an Epetra_Map object " 496 "that partitions elements of an index set " 497 "between processors. However, the union of the " 498 "index sets on different processors does not " 499 "contain all indices exactly once: the sum of " 500 "the number of entries the various processors " 501 "want to store locally is " 503 " whereas the total size of the object to be " 506 ". In other words, there are " 507 "either indices that are not spoken for " 508 "by any processor, or there are indices that are " 509 "claimed by multiple processors."));
514 return Epetra_Map (TrilinosWrappers::types::int_type(
size()),
515 TrilinosWrappers::types::int_type(
n_elements()),
517 #ifdef DEAL_II_WITH_MPI
518 Epetra_MpiComm(communicator));
520 Epetra_SerialComm());
524 std::vector<size_type> indices;
527 return Epetra_Map (TrilinosWrappers::types::int_type(-1),
528 TrilinosWrappers::types::int_type(
n_elements()),
531 reinterpret_cast<TrilinosWrappers::types::int_type *>(&indices[0])
535 #ifdef DEAL_II_WITH_MPI
536 Epetra_MpiComm(communicator));
538 Epetra_SerialComm());
558 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
void block_read(std::istream &in)
::ExceptionBase & ExcMessage(std::string arg1)
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
Epetra_Map make_trilinos_map(const MPI_Comm &communicator=MPI_COMM_WORLD, const bool overlapping=false) const
void block_write(std::ostream &out) const
std::vector< Range > ranges
T sum(const T &t, const MPI_Comm &mpi_communicator)
void subtract_set(const IndexSet &other)
#define Assert(cond, exc)
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
bool is_contiguous() const
void write(std::ostream &out) const
types::global_dof_index size_type
ElementIterator begin() const
void read(std::istream &in)
size_type n_elements() const