15#ifndef dealii_mpi_mpi_remote_point_evaluation_h
16#define dealii_mpi_mpi_remote_point_evaluation_h
30 template <
int dim,
int spacedim>
35 template <
int dim,
int spacedim>
54 template <
int dim,
int spacedim = dim>
128 DEAL_II_DEPRECATED_EARLY_WITH_COMMENT(
129 "Use the constructor with AdditionalData struct.")
131 const
double tolerance,
132 const
bool enforce_unique_mapping = false,
133 const
unsigned int rtree_level = 0,
134 const
std::function<
std::vector<
bool>()> &marked_vertices = {});
182 DistributedComputePointLocationsInternal<dim, spacedim> &data,
243 template <
typename T>
251 std::vector<std::pair<int, int>>
cells;
299 template <
typename T,
unsigned int n_components = 1>
302 std::vector<T> &output,
303 std::vector<T> &buffer,
304 const std::function<
void(
const ArrayView<T> &,
const CellData &)>
305 &evaluation_function,
306 const bool sort_data =
true)
const;
312 template <
typename T,
unsigned int n_components = 1>
315 const std::function<
void(
const ArrayView<T> &,
const CellData &)>
316 &evaluation_function,
317 const bool sort_data =
true)
const;
333 template <
typename T,
unsigned int n_components = 1>
336 const std::vector<T> &input,
337 std::vector<T> &buffer,
339 &evaluation_function,
340 const bool sort_data =
true)
const;
346 template <
typename T,
unsigned int n_components = 1>
349 const std::vector<T> &input,
351 &evaluation_function,
352 const bool sort_data =
true)
const;
358 const std::vector<unsigned int> &
407 const std::vector<unsigned int> &
416 const std::vector<unsigned int> &
537#ifdef DEAL_II_WITH_MPI
541 template <
typename T>
542 std::enable_if_t<Utilities::MPI::is_mpi_type<T> ==
false,
void>
544 const unsigned int rank,
545 const unsigned int tag,
547 std::vector<std::vector<char>> &buffers,
548 std::vector<MPI_Request> &requests)
550 requests.emplace_back(MPI_Request());
553 std::vector<T>(data.
data(), data.
data() + data.
size()),
false));
555 const int ierr = MPI_Isend(buffers.back().data(),
556 buffers.back().size(),
571 template <
typename T>
572 std::enable_if_t<Utilities::MPI::is_mpi_type<T> ==
true,
void>
574 const unsigned int rank,
575 const unsigned int tag,
577 std::vector<std::vector<char>> & ,
578 std::vector<MPI_Request> &requests)
580 requests.emplace_back(MPI_Request());
582 const int ierr = MPI_Isend(data.
data(),
598 template <
int rank_,
int dim,
typename T>
599 std::enable_if_t<Utilities::MPI::is_mpi_type<T> ==
true,
void>
601 const unsigned int rank,
602 const unsigned int tag,
604 std::vector<std::vector<char>> &buffers,
605 std::vector<MPI_Request> &requests)
617 template <
typename T>
618 std::enable_if_t<Utilities::MPI::is_mpi_type<T> ==
false,
void>
621 const MPI_Status &status,
622 std::vector<char> &buffer)
625 int ierr = MPI_Get_count(&status, MPI_CHAR, &message_length);
628 buffer.resize(message_length);
630 ierr = MPI_Recv(buffer.data(),
642 for (
unsigned int i = 0; i < data.
size(); ++i)
652 template <
typename T>
653 std::enable_if_t<Utilities::MPI::is_mpi_type<T> ==
true,
void>
656 const MPI_Status &status,
657 std::vector<char> & )
659 const auto ierr = MPI_Recv(data.
data(),
675 template <
int rank_,
int dim,
typename T>
676 std::enable_if_t<Utilities::MPI::is_mpi_type<T> ==
true,
void>
679 const MPI_Status &status,
680 std::vector<char> &buffer)
682 const ArrayView<T> data_(
reinterpret_cast<T *
>(data.data()),
692 template <
int dim,
int spacedim>
693 template <
typename T>
696 const unsigned int cell,
706 template <
int dim,
int spacedim>
707 template <
typename T,
unsigned int n_components>
710 std::vector<T> &output,
711 std::vector<T> &buffer,
713 &evaluation_function,
714 const bool sort_data)
const
716#ifndef DEAL_II_WITH_MPI
720 (void)evaluation_function;
726 const unsigned int my_rank =
730 output.resize(
point_ptrs.back() * n_components);
748 std::vector<MPI_Request> send_requests;
749 std::vector<std::vector<char>> send_buffers_packed;
750 std::vector<char> recv_buffer_packed;
753 evaluation_function(buffer_eval, *
cell_data);
758 const auto my_rank_local_recv_ptr =
761 if (my_rank_local_recv_ptr !=
recv_ranks.end())
763 const unsigned int my_rank_local_recv =
764 std::distance(
recv_ranks.begin(), my_rank_local_recv_ptr);
765 const unsigned int my_rank_local_send = std::distance(
768 const unsigned int start =
send_ptrs[my_rank_local_send];
769 const unsigned int end =
send_ptrs[my_rank_local_send + 1];
770 const unsigned int *recv_ptr =
776 if (start <= send_index && send_index < end)
778 for (
unsigned int c = 0; c < n_components; ++c)
779 output[recv_ptr[send_index - start] * n_components + c] =
780 buffer_eval[i * n_components + c];
783 for (
unsigned int c = 0; c < n_components; ++c)
784 buffer_send[send_index * n_components + c] =
785 buffer_eval[i * n_components + c];
791 for (
unsigned int c = 0; c < n_components; ++c)
793 buffer_eval[i * n_components + c];
798 send_buffers_packed.reserve(
send_ranks.size());
801 for (
unsigned int i = 0; i <
send_ranks.size(); ++i)
820 const auto my_rank_local_recv_ptr =
823 if (my_rank_local_recv_ptr !=
recv_ranks.end())
825 const unsigned int my_rank_local_recv =
826 std::distance(
recv_ranks.begin(), my_rank_local_recv_ptr);
827 const unsigned int my_rank_local_send = std::distance(
831 for (
unsigned int j =
recv_ptrs[my_rank_local_recv],
835 for (
unsigned int c = 0; c < n_components; ++c)
836 output[j * n_components + c] =
837 buffer_eval[k * n_components + c];
842 for (
unsigned int i = 0; i <
recv_ranks.size(); ++i)
849 int ierr = MPI_Probe(MPI_ANY_SOURCE,
860 const unsigned int j = std::distance(
recv_ranks.begin(), ptr);
866 (output.data() +
recv_ptrs[j] * n_components),
879 for (
unsigned int c = 0; c < n_components; ++c)
881 recv_buffer[k * n_components + c];
886 if (!send_requests.
empty())
888 const int ierr = MPI_Waitall(send_requests.
size(),
889 send_requests.
data(),
890 MPI_STATUSES_IGNORE);
897 template <
int dim,
int spacedim>
898 template <
typename T,
unsigned int n_components>
902 &evaluation_function,
903 const bool sort_data)
const
905 std::vector<T> output;
906 std::vector<T> buffer;
918 template <
int dim,
int spacedim>
919 template <
typename T,
unsigned int n_components>
922 const std::vector<T> &input,
923 std::vector<T> &buffer,
925 &evaluation_function,
926 const bool sort_data)
const
928#ifndef DEAL_II_WITH_MPI
932 (void)evaluation_function;
938 const unsigned int my_rank =
945 sort_data ? ((
point_ptrs.size() - 1) * n_components) :
959 const_cast<T *
>(input.data()),
963 std::vector<MPI_Request> send_requests;
964 std::vector<std::vector<char>> send_buffers_packed;
965 std::vector<char> recv_buffer_packed;
971 const auto my_rank_local_recv_ptr =
974 if (my_rank_local_recv_ptr !=
recv_ranks.end())
978 const unsigned int my_rank_local_recv =
979 std::distance(
recv_ranks.begin(), my_rank_local_recv_ptr);
980 const unsigned int my_rank_local_send = std::distance(
984 const unsigned int start =
recv_ptrs[my_rank_local_recv];
985 const unsigned int end =
recv_ptrs[my_rank_local_recv + 1];
986 const unsigned int *send_ptr =
988 for (
unsigned int i = 0, k = 0; i <
point_ptrs.size() - 1; ++i)
991 for (
unsigned int j =
point_ptrs[i]; j < next; ++j, ++k)
996 if (start <= recv_index && recv_index < end)
997 for (
unsigned int c = 0; c < n_components; ++c)
998 buffer_eval[send_ptr[recv_index - start] *
1000 c] = input[i * n_components + c];
1002 for (
unsigned int c = 0; c < n_components; ++c)
1003 buffer_send[recv_index * n_components + c] =
1004 input[i * n_components + c];
1010 for (
unsigned int i = 0, k = 0; i <
point_ptrs.size() - 1; ++i)
1013 for (
unsigned int c = 0; c < n_components; ++c)
1015 input[i * n_components + c];
1020 send_buffers_packed.reserve(
recv_ranks.size());
1023 for (
unsigned int i = 0; i <
recv_ranks.size(); ++i)
1035 send_buffers_packed,
1042 const auto my_rank_local_recv_ptr =
1045 if (my_rank_local_recv_ptr !=
recv_ranks.end())
1047 const unsigned int my_rank_local_recv =
1048 std::distance(
recv_ranks.begin(), my_rank_local_recv_ptr);
1049 const unsigned int my_rank_local_send = std::distance(
1053 for (
unsigned int j =
recv_ptrs[my_rank_local_recv],
1057 for (
unsigned int c = 0; c < n_components; ++c)
1058 buffer_eval[k * n_components + c] =
1059 input[j * n_components + c];
1063 for (
unsigned int i = 0; i <
send_ranks.size(); ++i)
1070 int ierr = MPI_Probe(MPI_ANY_SOURCE,
1082 const unsigned int j = std::distance(
send_ranks.begin(), ptr);
1088 (buffer.data() +
send_ptrs[j] * n_components),
1094 recv_buffer_packed);
1100 for (
unsigned int c = 0; c < n_components; ++c)
1102 recv_buffer[k * n_components + c];
1106 if (!send_requests.
empty())
1108 const int ierr = MPI_Waitall(send_requests.
size(),
1109 send_requests.
data(),
1110 MPI_STATUSES_IGNORE);
1115 evaluation_function(buffer_eval, *
cell_data);
1121 template <
int dim,
int spacedim>
1122 template <
typename T,
unsigned int n_components>
1125 const std::vector<T> &input,
1127 &evaluation_function,
1128 const bool sort_data)
const
1130 std::vector<T> buffer;
1133 evaluation_function,
value_type * data() const noexcept
Abstract base class for mapping classes.
virtual MPI_Comm get_communicator() const
ArrayView< T > get_data_view(const unsigned int cell, const ArrayView< T > &values) const
std::vector< std::pair< int, int > > cells
std_cxx20::ranges::iota_view< unsigned int, unsigned int > cell_indices() const
ArrayView< const Point< dim > > get_unit_points(const unsigned int cell) const
Triangulation< dim, spacedim >::active_cell_iterator get_active_cell_iterator(const unsigned int cell) const
const Triangulation< dim, spacedim > & triangulation
CellData(const Triangulation< dim, spacedim > &triangulation)
std::vector< unsigned int > reference_point_ptrs
std::vector< Point< dim > > reference_point_values
const std::vector< unsigned int > & get_send_permutation() const
std::vector< unsigned int > send_ptrs
std::unique_ptr< CellData > cell_data
bool all_points_found() const
unsigned int buffer_size_without_sorting
boost::signals2::connection tria_signal
std::vector< unsigned int > point_ptrs
RemotePointEvaluation(const AdditionalData &additional_data=AdditionalData())
std::vector< unsigned int > recv_permutation
const Triangulation< dim, spacedim > & get_triangulation() const
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 bool sort_data=true) const
SmartPointer< const Mapping< dim, spacedim > > mapping
const std::vector< unsigned int > & get_point_ptrs() const
std::vector< T > evaluate_and_process(const std::function< void(const ArrayView< T > &, const CellData &)> &evaluation_function, const bool sort_data=true) const
const std::vector< unsigned int > & get_inverse_recv_permutation() const
std::vector< unsigned int > send_permutation_inv
const AdditionalData additional_data
bool is_map_unique() const
std::vector< unsigned int > recv_ranks
std::vector< unsigned int > recv_permutation_inv
std::vector< unsigned int > recv_ptrs
void evaluate_and_process(std::vector< T > &output, std::vector< T > &buffer, const std::function< void(const ArrayView< T > &, const CellData &)> &evaluation_function, const bool sort_data=true) const
const Mapping< dim, spacedim > & get_mapping() const
const CellData & get_cell_data() const
std::vector< unsigned int > send_ranks
void process_and_evaluate(const std::vector< T > &input, const std::function< void(const ArrayView< const T > &, const CellData &)> &evaluation_function, const bool sort_data=true) const
bool all_points_found_flag
unsigned int buffer_size_with_sorting
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)
std::enable_if_t< Utilities::MPI::is_mpi_type< T >==false, void > recv_and_unpack(const ArrayView< T > &data, const MPI_Comm comm, const MPI_Status &status, std::vector< char > &buffer)
std::enable_if_t< Utilities::MPI::is_mpi_type< T >==false, void > pack_and_isend(const ArrayView< const T > &data, const unsigned int rank, const unsigned int tag, const MPI_Comm comm, std::vector< std::vector< char > > &buffers, std::vector< MPI_Request > &requests)
unsigned int this_mpi_process(const MPI_Comm mpi_communicator)
const MPI_Datatype mpi_type_id_for_type
size_t pack(const T &object, std::vector< char > &dest_buffer, const bool allow_compression=true)
T unpack(const std::vector< char > &buffer, const bool allow_compression=true)
constexpr T pow(const T base, const int iexp)
boost::integer_range< IncrementableType > iota_view
*braid_SplitCommworld & comm
AdditionalData(const double tolerance=1e-6, const bool enforce_unique_mapping=false, const unsigned int rtree_level=0, const std::function< std::vector< bool >()> &marked_vertices={})
std::function< std::vector< bool >()> marked_vertices
bool enforce_unique_mapping