42 #ifndef TPETRA_DETAILS_ALLREDUCEVIEW_HPP
43 #define TPETRA_DETAILS_ALLREDUCEVIEW_HPP
47 #include "Kokkos_Core.hpp"
48 #include "Teuchos_CommHelpers.hpp"
50 #include <type_traits>
58 template<
class ViewType>
59 struct view_is_cuda_uvm {
60 static constexpr
bool value =
61 #ifdef KOKKOS_ENABLE_CUDA
62 std::is_same<
typename ViewType::memory_space,
63 Kokkos::CudaUVMSpace>::value;
66 #endif // KOKKOS_ENABLE_CUDA
69 template<
class ViewType>
70 struct MakeContiguousBuffer {
71 static constexpr
bool is_contiguous_layout =
73 typename ViewType::array_layout,
74 Kokkos::LayoutLeft>::value ||
76 typename ViewType::array_layout,
77 Kokkos::LayoutRight>::value;
78 using contiguous_array_layout =
79 typename std::conditional<is_contiguous_layout,
80 typename ViewType::array_layout,
81 Kokkos::LayoutLeft>::type;
82 using contiguous_device_type =
83 typename std::conditional<
85 typename ViewType::memory_space,
86 Kokkos::HostSpace>::value,
87 typename ViewType::device_type,
88 Kokkos::HostSpace::device_type>::type;
89 using contiguous_buffer_type =
90 Kokkos::View<
typename ViewType::non_const_data_type,
91 contiguous_array_layout,
92 contiguous_device_type>;
94 static contiguous_array_layout
95 makeLayout (
const ViewType& view)
99 return contiguous_array_layout (view.extent (0), view.extent (1),
100 view.extent (2), view.extent (3),
101 view.extent (4), view.extent (5),
102 view.extent (6), view.extent (7));
105 static contiguous_buffer_type
106 make (
const ViewType& view)
108 using Kokkos::view_alloc;
109 using Kokkos::WithoutInitializing;
110 return contiguous_buffer_type
111 (view_alloc (view.label (), WithoutInitializing),
116 template<
class ViewType>
117 typename MakeContiguousBuffer<ViewType>::contiguous_buffer_type
118 makeContiguousBuffer (
const ViewType& view)
120 return MakeContiguousBuffer<ViewType>::make (view);
123 template<
class ValueType>
125 allReduceRawContiguous (ValueType output[],
126 const ValueType input[],
128 const Teuchos::Comm<int>& comm)
130 using Teuchos::outArg;
131 using Teuchos::REDUCE_SUM;
132 using Teuchos::reduceAll;
133 constexpr
size_t max_int = size_t (std::numeric_limits<int>::max ());
134 TEUCHOS_ASSERT( count <=
size_t (max_int) );
135 reduceAll<int, ValueType> (comm, REDUCE_SUM,
static_cast<int> (count),
147 template<
class InputViewType,
class OutputViewType>
150 const InputViewType& input,
151 const Teuchos::Comm<int>& comm)
175 const bool viewsAlias = output.data () == input.data ();
176 if (comm.getSize () == 1) {
187 const bool assumeMpiCanAccessBuffers =
188 ! view_is_cuda_uvm<OutputViewType>::value &&
189 ! view_is_cuda_uvm<InputViewType>::value &&
192 const bool needContiguousTemporaryBuffers =
193 ! assumeMpiCanAccessBuffers ||
194 ::Tpetra::Details::isInterComm (comm) ||
195 output.span_is_contiguous () || input.span_is_contiguous ();
196 if (needContiguousTemporaryBuffers) {
197 auto output_tmp = makeContiguousBuffer (output);
198 auto input_tmp = makeContiguousBuffer (input);
205 allReduceRawContiguous (output_tmp.data (), input_tmp.data (),
206 output_tmp.span (), comm);
210 allReduceRawContiguous (output.data (), input.data (),
211 output.span (), comm);
218 #endif // TPETRA_DETAILS_ALLREDUCEVIEW_HPP
void deep_copy(MultiVector< DS, DL, DG, DN > &dst, const MultiVector< SS, SL, SG, SN > &src)
Copy the contents of the MultiVector src into dst.
static void allReduceView(const OutputViewType &output, const InputViewType &input, const Teuchos::Comm< int > &comm)
All-reduce from input Kokkos::View to output Kokkos::View.
static bool assumeMpiIsCudaAware()
Whether to assume that MPI is CUDA aware.
Declaration of Tpetra::Details::isInterComm.
Declaration of Tpetra::Details::Behavior, a class that describes Tpetra's behavior.