16#ifndef dealii_mpi_mpi_remote_point_evaluation_h
17#define dealii_mpi_mpi_remote_point_evaluation_h
40 template <
int dim,
int spacedim = dim>
94 std::vector<std::pair<int, int>>
cells;
122 template <
typename T>
125 std::vector<T> &output,
126 std::vector<T> &buffer,
128 &evaluation_function)
const;
138 template <
typename T>
141 const std::vector<T> &input,
142 std::vector<T> & buffer,
144 &evaluation_function)
const;
150 const std::vector<unsigned int> &
297 template <
int dim,
int spacedim>
298 template <
typename T>
301 std::vector<T> &output,
302 std::vector<T> &buffer,
304 &evaluation_function)
const
306#ifndef DEAL_II_WITH_MPI
310 (void)evaluation_function;
315 output.resize(point_ptrs.back());
316 buffer.resize(send_permutation.size() * 2);
317 ArrayView<T> buffer_1(buffer.data(), buffer.size() / 2);
318 ArrayView<T> buffer_2(buffer.data() + buffer.size() / 2,
322 evaluation_function(buffer_1, cell_data);
325 for (
unsigned int i = 0; i < send_permutation.size(); ++i)
326 buffer_2[send_permutation[i]] = buffer_1[i];
329 std::map<unsigned int, std::vector<char>> temp_map;
331 std::vector<MPI_Request> requests;
332 requests.reserve(send_ranks.size());
334 const unsigned int my_rank =
337 std::map<unsigned int, std::vector<T>> temp_recv_map;
339 for (
unsigned int i = 0; i < send_ranks.size(); ++i)
341 if (send_ranks[i] == my_rank)
344 temp_recv_map[my_rank] =
345 std::vector<T>(buffer_2.
begin() + send_ptrs[i],
346 buffer_2.
begin() + send_ptrs[i + 1]);
350 temp_map[send_ranks[i]] =
352 buffer_2.
begin() + send_ptrs[i + 1]),
355 auto &buffer = temp_map[send_ranks[i]];
357 requests.push_back(MPI_Request());
359 const int ierr = MPI_Isend(buffer.data(),
369 for (
const auto recv_rank : recv_ranks)
371 if (recv_rank == my_rank)
375 int ierr = MPI_Probe(MPI_ANY_SOURCE,
382 ierr = MPI_Get_count(&status, MPI_CHAR, &message_length);
385 std::vector<char> buffer(message_length);
387 ierr = MPI_Recv(buffer.data(),
396 temp_recv_map[status.MPI_SOURCE] =
397 Utilities::unpack<std::vector<T>>(buffer,
false);
402 MPI_Waitall(requests.size(), requests.data(), MPI_STATUSES_IGNORE);
406 auto it = recv_permutation.begin();
407 for (
const auto &j : temp_recv_map)
408 for (
const auto &i : j.second)
417 template <
int dim,
int spacedim>
418 template <
typename T>
421 const std::vector<T> &input,
422 std::vector<T> & buffer,
424 &evaluation_function)
const
426#ifndef DEAL_II_WITH_MPI
430 (void)evaluation_function;
435 const auto &ptr = this->get_point_ptrs();
437 std::map<unsigned int, std::vector<T>> temp_recv_map;
439 for (
unsigned int i = 0; i < recv_ranks.size(); ++i)
440 temp_recv_map[recv_ranks[i]].resize(recv_ptrs[i + 1] - recv_ptrs[i]);
442 const unsigned int my_rank =
449 for (
auto &j : temp_recv_map)
450 i += j.second.size();
458 std::vector<T> buffer_(ptr.back());
459 for (
unsigned int i = 0, c = 0; i < ptr.size() - 1; ++i)
461 const auto n_entries = ptr[i + 1] - ptr[i];
463 for (
unsigned int j = 0; j < n_entries; ++j, ++c)
464 buffer_[c] = input[i];
468 auto it = recv_permutation.begin();
469 for (
auto &j : temp_recv_map)
470 for (
auto &i : j.second)
478 buffer.resize(send_permutation.size() * 2);
479 ArrayView<T> buffer_1(buffer.data(), buffer.size() / 2);
480 ArrayView<T> buffer_2(buffer.data() + buffer.size() / 2,
484 std::map<unsigned int, std::vector<char>> temp_map;
486 std::vector<MPI_Request> requests;
487 requests.reserve(recv_ranks.size());
489 for (
const auto recv_rank : recv_ranks)
491 if (recv_rank == my_rank)
494 temp_map[recv_rank] =
497 auto &buffer_send = temp_map[recv_rank];
499 requests.push_back(MPI_Request());
501 const int ierr = MPI_Isend(buffer_send.data(),
511 for (
unsigned int i = 0; i < send_ranks.size(); ++i)
513 if (send_ranks[i] == my_rank)
515 const auto &buffer_send = temp_recv_map[send_ranks[i]];
517 const unsigned int j = std::distance(send_ranks.begin(),
518 std::find(send_ranks.begin(),
523 send_ptrs[j + 1] - send_ptrs[j]);
525 for (
unsigned int i = send_ptrs[j], c = 0; i < send_ptrs[j + 1];
527 buffer_1[i] = buffer_send[c];
533 int ierr = MPI_Probe(MPI_ANY_SOURCE,
540 ierr = MPI_Get_count(&status, MPI_CHAR, &message_length);
543 std::vector<char> recv_buffer(message_length);
545 ierr = MPI_Recv(recv_buffer.data(),
554 const auto recv_buffer_unpacked =
555 Utilities::unpack<std::vector<T>>(recv_buffer,
false);
558 std::find(send_ranks.begin(), send_ranks.end(), status.MPI_SOURCE);
562 const unsigned int j = std::distance(send_ranks.begin(), ptr);
565 send_ptrs[j + 1] - send_ptrs[j]);
567 for (
unsigned int i = send_ptrs[j], c = 0; i < send_ptrs[j + 1];
572 buffer_1[i] = recv_buffer_unpacked[c];
577 MPI_Waitall(requests.size(), requests.data(), MPI_STATUSES_IGNORE);
581 for (
unsigned int i = 0; i < send_permutation.size(); ++i)
582 buffer_2[i] = buffer_1[send_permutation[i]];
585 evaluation_function(buffer_2, cell_data);
Abstract base class for mapping classes.
virtual MPI_Comm get_communicator() const
void evaluate_and_process(std::vector< T > &output, std::vector< T > &buffer, const std::function< void(const ArrayView< T > &, const CellData &)> &evaluation_function) const
std::vector< unsigned int > send_ptrs
bool all_points_found() const
const bool enforce_unique_mapping
boost::signals2::connection tria_signal
std::vector< unsigned int > point_ptrs
std::vector< unsigned int > recv_permutation
const Triangulation< dim, spacedim > & get_triangulation() const
SmartPointer< const Mapping< dim, spacedim > > mapping
void process_and_evaluate(const std::vector< T > &input, std::vector< T > &buffer, const std::function< void(const ArrayView< const T > &, const CellData &)> &evaluation_function) const
const std::vector< unsigned int > & get_point_ptrs() const
const unsigned int rtree_level
bool is_map_unique() const
std::vector< unsigned int > recv_ranks
std::vector< unsigned int > recv_ptrs
const Mapping< dim, spacedim > & get_mapping() const
const std::function< std::vector< bool >()> marked_vertices
std::vector< unsigned int > send_ranks
bool all_points_found_flag
void reinit(const std::vector< Point< spacedim > > &points, const Triangulation< dim, spacedim > &tria, const Mapping< dim, spacedim > &mapping)
SmartPointer< const Triangulation< dim, spacedim > > tria
bool point_found(const unsigned int i) const
std::vector< unsigned int > send_permutation
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcNotImplemented()
static ::ExceptionBase & ExcNeedsMPI()
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
#define AssertThrowMPI(error_code)
#define AssertIndexRange(index, range)
unsigned int this_mpi_process(const MPI_Comm &mpi_communicator)
size_t pack(const T &object, std::vector< char > &dest_buffer, const bool allow_compression=true)
std::vector< std::pair< int, int > > cells
std::vector< unsigned int > reference_point_ptrs
std::vector< Point< dim > > reference_point_values
const ::Triangulation< dim, spacedim > & tria