16 #include <deal.II/base/partitioner.h> 18 DEAL_II_NAMESPACE_OPEN
27 local_range_data (
std::pair<
types::global_dof_index,
types::global_dof_index> (0, 0)),
28 n_ghost_indices_data (0),
29 n_import_indices_data (0),
32 communicator (MPI_COMM_SELF),
33 have_ghost_indices (false)
42 locally_owned_range_data (size),
43 local_range_data (
std::pair<
types::global_dof_index,
types::global_dof_index> (0, size)),
44 n_ghost_indices_data (0),
45 n_import_indices_data (0),
48 communicator (MPI_COMM_SELF),
49 have_ghost_indices (false)
60 const MPI_Comm communicator_in)
62 global_size (static_cast<
types::global_dof_index>(locally_owned_indices.size())),
63 n_ghost_indices_data (0),
64 n_import_indices_data (0),
67 communicator (communicator_in),
68 have_ghost_indices (false)
77 const MPI_Comm communicator_in)
79 global_size (static_cast<
types::global_dof_index>(locally_owned_indices.size())),
80 n_ghost_indices_data (0),
81 n_import_indices_data (0),
84 communicator (communicator_in),
85 have_ghost_indices (false)
94 const IndexSet &read_write_vector_index_set,
95 const MPI_Comm &communicator_in)
121 ExcMessage (
"The index set specified in locally_owned_indices " 122 "is not contiguous."));
125 local_range_data = std::pair<types::global_dof_index, types::global_dof_index>
130 static_cast<types::global_dof_index>(std::numeric_limits<unsigned int>::max()),
131 ExcMessage(
"Index overflow: This class supports at most 2^32-1 locally owned vector entries"));
159 ExcMessage(
"Index overflow: This class supports at most 2^32-1 ghost elements"));
173 #ifdef DEAL_II_WITH_MPI 182 std::vector<types::global_dof_index> first_index (
n_procs+1);
186 int ierr = MPI_Bcast(&first_index[0], 1, DEAL_II_DOF_INDEX_MPI_TYPE,
192 DEAL_II_DOF_INDEX_MPI_TYPE, &first_index[1], 1,
201 unsigned int first_proc_with_nonzero_dofs = 0;
202 for (
unsigned int i=0; i<
n_procs; ++i)
203 if (first_index[i+1]>0)
205 first_proc_with_nonzero_dofs = i;
208 for (
unsigned int i=first_proc_with_nonzero_dofs+1; i<
n_procs; ++i)
209 if (first_index[i] == 0)
210 first_index[i] = first_index[i-1];
223 unsigned int n_ghost_targets = 0;
231 unsigned int current_proc = 0;
234 while (current_index >= first_index[current_proc+1])
236 std::vector<std::pair<unsigned int, unsigned int> > ghost_targets_temp
237 (1, std::pair<unsigned int, unsigned int>(current_proc, 0));
242 current_index = expanded_ghost_indices[iterator];
243 while (current_index >= first_index[current_proc+1])
246 if ( ghost_targets_temp[n_ghost_targets-1].first < current_proc)
248 ghost_targets_temp[n_ghost_targets-1].second =
249 iterator - ghost_targets_temp[n_ghost_targets-1].second;
250 ghost_targets_temp.push_back(std::pair<
unsigned int,
251 unsigned int>(current_proc,iterator));
255 ghost_targets_temp[n_ghost_targets-1].second =
261 std::vector<int> send_buffer (
n_procs, 0);
262 std::vector<int> receive_buffer (
n_procs, 0);
263 for (
unsigned int i=0; i<n_ghost_targets; i++)
266 const int ierr = MPI_Alltoall (&send_buffer[0], 1, MPI_INT, &receive_buffer[0], 1,
271 std::vector<std::pair<unsigned int,unsigned int> > import_targets_temp;
273 for (
unsigned int i=0; i<
n_procs; i++)
274 if (receive_buffer[i] > 0)
277 import_targets_temp.push_back(std::pair<
unsigned int,
278 unsigned int> (i, receive_buffer[i]));
287 unsigned int current_index_start = 0;
291 const int ierr = MPI_Irecv (&expanded_import_indices[current_index_start],
293 DEAL_II_DOF_INDEX_MPI_TYPE,
303 current_index_start = 0;
304 for (
unsigned int i=0; i<n_ghost_targets; i++)
306 const int ierr = MPI_Send (&expanded_ghost_indices[current_index_start],
315 if (import_requests.size()>0)
317 const int ierr = MPI_Waitall (import_requests.size(),
319 MPI_STATUSES_IGNORE);
327 std::vector<std::pair<unsigned int,unsigned int> > compressed_import_indices;
338 if (new_index == last_index+1)
339 compressed_import_indices.back().second++;
342 compressed_import_indices.push_back
343 (std::pair<unsigned int,unsigned int>(new_index,new_index+1));
345 last_index = new_index;
372 #ifdef DEAL_II_WITH_MPI 375 int communicators_same = 0;
377 &communicators_same);
379 if (!(communicators_same == MPI_IDENT ||
380 communicators_same == MPI_CONGRUENT))
418 DEAL_II_NAMESPACE_CLOSE
std::vector< std::pair< unsigned int, unsigned int > > import_indices_data
bool is_globally_compatible(const Partitioner &part) const
static const unsigned int invalid_unsigned_int
#define AssertDimension(dim1, dim2)
virtual void reinit(const IndexSet &vector_space_vector_index_set, const IndexSet &read_write_vector_index_set, const MPI_Comm &communicator)
bool is_compatible(const Partitioner &part) const
size_type nth_index_in_set(const unsigned int local_index) const
IndexSet ghost_indices_data
types::global_dof_index global_size
#define AssertIndexRange(index, range)
unsigned int n_ghost_indices_data
#define AssertThrow(cond, exc)
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
void set_owned_indices(const IndexSet &locally_owned_indices)
static ::ExceptionBase & ExcMessage(std::string arg1)
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)
IndexSet locally_owned_range_data
bool job_supports_mpi() 1
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)
types::global_dof_index size() const
std::size_t memory_consumption() const
#define AssertThrowMPI(error_code)
bool is_contiguous() const
T min(const T &t, const MPI_Comm &mpi_communicator)
std::vector< std::pair< unsigned int, unsigned int > > import_targets_data
std::vector< std::pair< unsigned int, unsigned int > > ghost_targets_data
unsigned int this_mpi_process(const MPI_Comm &mpi_communicator)
static ::ExceptionBase & ExcNotImplemented()
std::pair< types::global_dof_index, types::global_dof_index > local_range_data
const types::global_dof_index invalid_dof_index
unsigned int n_import_indices_data
size_type n_elements() const
void set_ghost_indices(const IndexSet &ghost_indices)
static ::ExceptionBase & ExcInternalError()