16#ifndef CABANA_GRID_HALO_HPP
17#define CABANA_GRID_HALO_HPP
24#include <Kokkos_Core.hpp>
25#include <Kokkos_Profiling_ScopedRegion.hpp>
43template <std::
size_t NumSpaceDim>
58 setNeighbors(
const std::vector<std::array<int, num_space_dim>>& neighbors )
60 _neighbors = neighbors;
70 std::vector<std::array<int, num_space_dim>> _neighbors;
75template <std::
size_t NumSpaceDim>
81class NodeHaloPattern<3> :
public HaloPattern<3>
87 std::vector<std::array<int, 3>> neighbors;
88 neighbors.reserve( 26 );
89 for (
int i = -1; i < 2; ++i )
90 for (
int j = -1; j < 2; ++j )
91 for (
int k = -1; k < 2; ++k )
92 if ( !( i == 0 && j == 0 && k == 0 ) )
93 neighbors.push_back( { i, j, k } );
101class NodeHaloPattern<2> :
public HaloPattern<2>
107 std::vector<std::array<int, 2>> neighbors;
108 neighbors.reserve( 8 );
109 for (
int i = -1; i < 2; ++i )
110 for (
int j = -1; j < 2; ++j )
111 if ( !( i == 0 && j == 0 ) )
112 neighbors.push_back( { i, j } );
119template <std::
size_t NumSpaceDim>
125class FaceHaloPattern<3> :
public HaloPattern<3>
131 std::vector<std::array<int, 3>> neighbors = {
132 { -1, 0, 0 }, { 1, 0, 0 }, { 0, -1, 0 },
133 { 0, 1, 0 }, { 0, 0, -1 }, { 0, 0, 1 } };
141class FaceHaloPattern<2> :
public HaloPattern<2>
147 std::vector<std::array<int, 2>> neighbors = {
148 { -1, 0 }, { 1, 0 }, { 0, -1 }, { 0, 1 } };
156namespace ScatterReduce
199template <
class MemorySpace>
215 template <
class Pattern,
class... ArrayTypes>
216 Halo(
const Pattern& pattern,
const int width,
const ArrayTypes&... arrays )
219 const std::size_t num_space_dim = Pattern::num_space_dim;
225 auto neighbor_id = [](
const std::array<int, num_space_dim>& ijk )
228 for ( std::size_t d = 1; d < num_space_dim; ++d )
229 id += num_space_dim *
id + ijk[d];
235 auto flip_id = [=](
const std::array<int, num_space_dim>& ijk )
237 std::array<int, num_space_dim> flip_ijk;
238 for ( std::size_t d = 0; d < num_space_dim; ++d )
239 flip_ijk[d] = -ijk[d];
246 auto neighbors = pattern.getNeighbors();
247 for (
const auto& n : neighbors )
250 int rank = local_grid->neighborRank( n );
256 _neighbor_ranks.push_back( rank );
260 _send_tags.push_back( neighbor_id( n ) );
264 _receive_tags.push_back( neighbor_id( flip_id( n ) ) );
272 _ghosted_steering, arrays... );
287 template <
class ExecutionSpace,
class... ArrayTypes>
288 void gather(
const ExecutionSpace& exec_space,
289 const ArrayTypes&... arrays )
const
291 Kokkos::Profiling::ScopedRegion region(
"Cabana::Grid::gather" );
294 int num_n = _neighbor_ranks.size();
299 auto comm =
getComm( arrays... );
302 std::vector<MPI_Request> requests( 2 * num_n, MPI_REQUEST_NULL );
306 const int mpi_tag = 1234;
309 for (
int n = 0; n < num_n; ++n )
312 if ( 0 < _ghosted_buffers[n].
size() )
314 MPI_Irecv( _ghosted_buffers[n].data(),
315 _ghosted_buffers[n].
size(), MPI_BYTE,
316 _neighbor_ranks[n], mpi_tag + _receive_tags[n], comm,
322 for (
int n = 0; n < num_n; ++n )
325 if ( 0 < _owned_buffers[n].
size() )
328 packBuffer( exec_space, _owned_buffers[n], _owned_steering[n],
332 MPI_Isend( _owned_buffers[n].data(), _owned_buffers[n].
size(),
333 MPI_BYTE, _neighbor_ranks[n],
334 mpi_tag + _send_tags[n], comm,
335 &requests[num_n + n] );
340 bool unpack_complete =
false;
341 while ( !unpack_complete )
344 int unpack_index = MPI_UNDEFINED;
345 MPI_Waitany( num_n, requests.data(), &unpack_index,
349 if ( MPI_UNDEFINED == unpack_index )
351 unpack_complete =
true;
358 _ghosted_buffers[unpack_index],
359 _ghosted_steering[unpack_index],
365 MPI_Waitall( num_n, requests.data() + num_n, MPI_STATUSES_IGNORE );
375 template <
class ExecutionSpace,
class ReduceOp,
class... ArrayTypes>
376 void scatter(
const ExecutionSpace& exec_space,
const ReduceOp& reduce_op,
377 const ArrayTypes&... arrays )
const
379 Kokkos::Profiling::ScopedRegion region(
"Cabana::Grid::scatter" );
382 int num_n = _neighbor_ranks.size();
387 auto comm =
getComm( arrays... );
390 std::vector<MPI_Request> requests( 2 * num_n, MPI_REQUEST_NULL );
394 const int mpi_tag = 2345;
397 for (
int n = 0; n < num_n; ++n )
400 if ( 0 < _owned_buffers[n].
size() )
402 MPI_Irecv( _owned_buffers[n].data(), _owned_buffers[n].
size(),
403 MPI_BYTE, _neighbor_ranks[n],
404 mpi_tag + _receive_tags[n], comm, &requests[n] );
409 for (
int n = 0; n < num_n; ++n )
412 if ( 0 < _ghosted_buffers[n].
size() )
416 _ghosted_steering[n], arrays.view()... );
419 MPI_Isend( _ghosted_buffers[n].data(),
420 _ghosted_buffers[n].
size(), MPI_BYTE,
421 _neighbor_ranks[n], mpi_tag + _send_tags[n], comm,
422 &requests[num_n + n] );
427 bool unpack_complete =
false;
428 while ( !unpack_complete )
431 int unpack_index = MPI_UNDEFINED;
432 MPI_Waitany( num_n, requests.data(), &unpack_index,
436 if ( MPI_UNDEFINED == unpack_index )
438 unpack_complete =
true;
445 _owned_buffers[unpack_index],
446 _owned_steering[unpack_index], arrays.view()... );
450 MPI_Waitall( num_n, requests.data() + num_n, MPI_STATUSES_IGNORE );
456 template <
class Array_t>
457 MPI_Comm
getComm(
const Array_t& array )
const
459 return array.layout()->localGrid()->globalGrid().comm();
463 template <
class Array_t,
class... ArrayTypes>
464 MPI_Comm
getComm(
const Array_t& array,
const ArrayTypes&... arrays )
const
471 MPI_Comm_compare( comm,
getComm( arrays... ), &result );
472 if ( result != MPI_IDENT && result != MPI_CONGRUENT )
473 throw std::runtime_error(
"Cabana::Grid::Halo::getComm: Arrays "
474 "have different communicators" );
481 template <
class Array_t>
484 return array.layout()->localGrid();
489 template <
class Array_t,
class... ArrayTypes>
496 if ( local_grid->haloCellWidth() !=
497 array.layout()->localGrid()->haloCellWidth() )
499 throw std::runtime_error(
"Cabana::Grid::Halo::getlocalGrid: "
500 "Arrays have different halo widths" );
507 template <
class DecompositionTag, std::size_t NumSpaceDim,
511 const std::array<int, NumSpaceDim>& nid,
512 std::vector<Kokkos::View<char*, memory_space>>& buffers,
513 std::vector<Kokkos::View<int**, memory_space>>& steering,
514 const ArrayTypes&... arrays )
517 const std::size_t num_array =
sizeof...( ArrayTypes );
520 std::array<std::size_t, num_array> value_byte_sizes = {
521 sizeof(
typename ArrayTypes::value_type )... };
525 std::array<IndexSpace<NumSpaceDim + 1>, num_array> spaces = {
526 ( arrays.layout()->sharedIndexSpace( decomposition_tag, nid,
531 int buffer_bytes = 0;
532 int buffer_num_element = 0;
533 for ( std::size_t a = 0; a < num_array; ++a )
535 buffer_bytes += value_byte_sizes[a] * spaces[a].size();
536 buffer_num_element += spaces[a].size();
542 Kokkos::View<char*, memory_space>(
"halo_buffer", buffer_bytes ) );
545 steering.push_back( Kokkos::View<int**, memory_space>(
546 "steering", buffer_num_element, 3 + NumSpaceDim ) );
550 buffer_num_element, steering );
554 template <std::
size_t NumArray>
557 const std::array<std::size_t, NumArray>& value_byte_sizes,
558 const int buffer_bytes,
const int buffer_num_element,
559 std::vector<Kokkos::View<int**, memory_space>>& steering )
565 Kokkos::create_mirror_view( Kokkos::HostSpace(), steering.back() );
566 int elem_counter = 0;
567 int byte_counter = 0;
568 for ( std::size_t a = 0; a < NumArray; ++a )
570 for (
int i = spaces[a].min( 0 ); i < spaces[a].max( 0 ); ++i )
572 for (
int j = spaces[a].min( 1 ); j < spaces[a].max( 1 ); ++j )
574 for (
int k = spaces[a].min( 2 ); k < spaces[a].max( 2 );
577 for (
int l = spaces[a].min( 3 );
578 l < spaces[a].max( 3 ); ++l )
581 host_steering( elem_counter, 0 ) = byte_counter;
584 host_steering( elem_counter, 1 ) = a;
587 host_steering( elem_counter, 2 ) = i;
588 host_steering( elem_counter, 3 ) = j;
589 host_steering( elem_counter, 4 ) = k;
590 host_steering( elem_counter, 5 ) = l;
596 byte_counter += value_byte_sizes[a];
604 if ( byte_counter != buffer_bytes )
605 throw std::logic_error(
"Cabana::Grid::Halo::buildSteeringVector: "
606 "Steering vector contains different number "
607 "of bytes than buffer" );
608 if ( elem_counter != buffer_num_element )
609 throw std::logic_error(
"Cabana::Grid::Halo::buildSteeringVector: "
610 "Steering vector contains different number "
611 "of elements than buffer" );
614 Kokkos::deep_copy( steering.back(), host_steering );
618 template <std::
size_t NumArray>
621 const std::array<std::size_t, NumArray>& value_byte_sizes,
622 const int buffer_bytes,
const int buffer_num_element,
623 std::vector<Kokkos::View<int**, memory_space>>& steering )
629 Kokkos::create_mirror_view( Kokkos::HostSpace(), steering.back() );
630 int elem_counter = 0;
631 int byte_counter = 0;
632 for ( std::size_t a = 0; a < NumArray; ++a )
634 for (
int i = spaces[a].min( 0 ); i < spaces[a].max( 0 ); ++i )
636 for (
int j = spaces[a].min( 1 ); j < spaces[a].max( 1 ); ++j )
638 for (
int l = spaces[a].min( 2 ); l < spaces[a].max( 2 );
642 host_steering( elem_counter, 0 ) = byte_counter;
645 host_steering( elem_counter, 1 ) = a;
648 host_steering( elem_counter, 2 ) = i;
649 host_steering( elem_counter, 3 ) = j;
650 host_steering( elem_counter, 4 ) = l;
656 byte_counter += value_byte_sizes[a];
663 if ( byte_counter != buffer_bytes )
664 throw std::logic_error(
"Cabana::Grid::Halo::buildSteeringVector: "
665 "Steering vector contains different number "
666 "of bytes than buffer" );
667 if ( elem_counter != buffer_num_element )
668 throw std::logic_error(
"Cabana::Grid::Halo::buildSteeringVector: "
669 "Steering vector contains different number "
670 "of elements than buffer" );
673 Kokkos::deep_copy( steering.back(), host_steering );
678 template <
class ArrayView>
679 KOKKOS_INLINE_FUNCTION
static std::enable_if_t<4 == ArrayView::rank, void>
681 const Kokkos::View<int**, memory_space>& steering,
682 const int element_idx,
const ArrayView& array_view )
684 const char* elem_ptr =
reinterpret_cast<const char*
>( &array_view(
685 steering( element_idx, 2 ), steering( element_idx, 3 ),
686 steering( element_idx, 4 ), steering( element_idx, 5 ) ) );
687 for ( std::size_t b = 0; b <
sizeof(
typename ArrayView::value_type );
690 buffer( steering( element_idx, 0 ) + b ) = *( elem_ptr + b );
696 template <
class ArrayView>
697 KOKKOS_INLINE_FUNCTION
static std::enable_if_t<3 == ArrayView::rank, void>
699 const Kokkos::View<int**, memory_space>& steering,
700 const int element_idx,
const ArrayView& array_view )
702 const char* elem_ptr =
reinterpret_cast<const char*
>(
703 &array_view( steering( element_idx, 2 ), steering( element_idx, 3 ),
704 steering( element_idx, 4 ) ) );
705 for ( std::size_t b = 0; b <
sizeof(
typename ArrayView::value_type );
708 buffer( steering( element_idx, 0 ) + b ) = *( elem_ptr + b );
713 template <
class... ArrayViews>
714 KOKKOS_INLINE_FUNCTION
static void
715 packArray(
const Kokkos::View<char*, memory_space>& buffer,
716 const Kokkos::View<int**, memory_space>& steering,
717 const int element_idx,
718 const std::integral_constant<std::size_t, 0>,
722 if ( 0 == steering( element_idx, 1 ) )
728 template <std::size_t N,
class... ArrayViews>
729 KOKKOS_INLINE_FUNCTION
static void
730 packArray(
const Kokkos::View<char*, memory_space>& buffer,
731 const Kokkos::View<int**, memory_space>& steering,
732 const int element_idx,
733 const std::integral_constant<std::size_t, N>,
737 if ( N == steering( element_idx, 1 ) )
745 packArray( buffer, steering, element_idx,
746 std::integral_constant<std::size_t, N - 1>(), array_views );
750 template <
class ExecutionSpace,
class... ArrayViews>
752 const Kokkos::View<char*, memory_space>& buffer,
753 const Kokkos::View<int**, memory_space>& steering,
754 ArrayViews... array_views )
const
757 Kokkos::parallel_for(
758 "Cabana::Grid::Halo::pack_buffer",
759 Kokkos::RangePolicy<ExecutionSpace>( exec_space, 0,
760 steering.extent( 0 ) ),
761 KOKKOS_LAMBDA(
const int i ) {
764 std::integral_constant<std::size_t,
765 sizeof...( ArrayViews ) - 1>(),
773 KOKKOS_INLINE_FUNCTION
static void
776 array_val += buffer_val;
781 KOKKOS_INLINE_FUNCTION
static void
784 if ( buffer_val < array_val )
785 array_val = buffer_val;
790 KOKKOS_INLINE_FUNCTION
static void
793 if ( buffer_val > array_val )
794 array_val = buffer_val;
799 KOKKOS_INLINE_FUNCTION
static void
802 array_val = buffer_val;
807 template <
class ReduceOp,
class ArrayView>
808 KOKKOS_INLINE_FUNCTION
static std::enable_if_t<4 == ArrayView::rank, void>
810 const Kokkos::View<char*, memory_space>& buffer,
811 const Kokkos::View<int**, memory_space>& steering,
812 const int element_idx,
const ArrayView& array_view )
814 typename ArrayView::value_type elem;
815 char* elem_ptr =
reinterpret_cast<char*
>( &elem );
816 for ( std::size_t b = 0; b <
sizeof(
typename ArrayView::value_type );
819 *( elem_ptr + b ) = buffer( steering( element_idx, 0 ) + b );
822 array_view( steering( element_idx, 2 ),
823 steering( element_idx, 3 ),
824 steering( element_idx, 4 ),
825 steering( element_idx, 5 ) ) );
830 template <
class ReduceOp,
class ArrayView>
831 KOKKOS_INLINE_FUNCTION
static std::enable_if_t<3 == ArrayView::rank, void>
833 const Kokkos::View<char*, memory_space>& buffer,
834 const Kokkos::View<int**, memory_space>& steering,
835 const int element_idx,
const ArrayView& array_view )
837 typename ArrayView::value_type elem;
838 char* elem_ptr =
reinterpret_cast<char*
>( &elem );
839 for ( std::size_t b = 0; b <
sizeof(
typename ArrayView::value_type );
842 *( elem_ptr + b ) = buffer( steering( element_idx, 0 ) + b );
845 array_view( steering( element_idx, 2 ),
846 steering( element_idx, 3 ),
847 steering( element_idx, 4 ) ) );
851 template <
class ReduceOp,
class... ArrayViews>
852 KOKKOS_INLINE_FUNCTION
static void
854 const Kokkos::View<char*, memory_space>& buffer,
855 const Kokkos::View<int**, memory_space>& steering,
856 const int element_idx,
857 const std::integral_constant<std::size_t, 0>,
861 if ( 0 == steering( element_idx, 1 ) )
867 template <
class ReduceOp, std::size_t N,
class... ArrayViews>
868 KOKKOS_INLINE_FUNCTION
static void
870 const Kokkos::View<char*, memory_space>& buffer,
871 const Kokkos::View<int**, memory_space>& steering,
872 const int element_idx,
873 const std::integral_constant<std::size_t, N>,
877 if ( N == steering( element_idx, 1 ) )
885 unpackArray( reduce_op, buffer, steering, element_idx,
886 std::integral_constant<std::size_t, N - 1>(),
891 template <
class ExecutionSpace,
class ReduceOp,
class... ArrayViews>
893 const ExecutionSpace& exec_space,
894 const Kokkos::View<char*, memory_space>& buffer,
895 const Kokkos::View<int**, memory_space>& steering,
896 ArrayViews... array_views )
const
899 Kokkos::parallel_for(
900 "Cabana::Grid::Halo::unpack_buffer",
901 Kokkos::RangePolicy<ExecutionSpace>( exec_space, 0,
902 steering.extent( 0 ) ),
903 KOKKOS_LAMBDA(
const int i ) {
905 reduce_op, buffer, steering, i,
906 std::integral_constant<std::size_t,
907 sizeof...( ArrayViews ) - 1>(),
914 std::vector<int> _neighbor_ranks;
917 std::vector<int> _send_tags;
920 std::vector<int> _receive_tags;
923 std::vector<Kokkos::View<char*, memory_space>> _owned_buffers;
926 std::vector<Kokkos::View<char*, memory_space>> _ghosted_buffers;
929 std::vector<Kokkos::View<int**, memory_space>> _owned_steering;
932 std::vector<Kokkos::View<int**, memory_space>> _ghosted_steering;
939template <
class ArrayT,
class... Types>
943 using type =
typename ArrayT::memory_space;
954template <
class Pattern,
class... ArrayTypes>
956 const ArrayTypes&... arrays )
959 return std::make_shared<Halo<memory_space>>( pattern, width, arrays... );
auto createHalo(const Pattern &pattern, const int width, const ArrayTypes &... arrays)
Halo creation function.
Definition Cabana_Grid_Halo.hpp:955
Pack variadic template parameters for device capture.
Definition Cabana_Grid_Halo.hpp:120
Base halo exchange pattern class.
Definition Cabana_Grid_Halo.hpp:45
std::vector< std::array< int, num_space_dim > > getNeighbors() const
Get the neighbors that are in the halo pattern.
Definition Cabana_Grid_Halo.hpp:64
void setNeighbors(const std::vector< std::array< int, num_space_dim > > &neighbors)
Assign the neighbors that are in the halo pattern.
Definition Cabana_Grid_Halo.hpp:58
static constexpr std::size_t num_space_dim
Spatial dimension.
Definition Cabana_Grid_Halo.hpp:48
static KOKKOS_INLINE_FUNCTION void unpackOp(ScatterReduce::Sum, const T &buffer_val, T &array_val)
Reduce an element into the buffer. Sum reduction.
Definition Cabana_Grid_Halo.hpp:774
MemorySpace memory_space
Memory space.
Definition Cabana_Grid_Halo.hpp:204
void scatter(const ExecutionSpace &exec_space, const ReduceOp &reduce_op, const ArrayTypes &... arrays) const
Scatter data from our ghosts to their owners using the given type of reduce operation.
Definition Cabana_Grid_Halo.hpp:376
static KOKKOS_INLINE_FUNCTION void packArray(const Kokkos::View< char *, memory_space > &buffer, const Kokkos::View< int **, memory_space > &steering, const int element_idx, const std::integral_constant< std::size_t, 0 >, const Cabana::ParameterPack< ArrayViews... > &array_views)
Pack an array into a buffer.
Definition Cabana_Grid_Halo.hpp:715
static KOKKOS_INLINE_FUNCTION std::enable_if_t< 4==ArrayView::rank, void > unpackElement(const ReduceOp &reduce_op, const Kokkos::View< char *, memory_space > &buffer, const Kokkos::View< int **, memory_space > &steering, const int element_idx, const ArrayView &array_view)
Definition Cabana_Grid_Halo.hpp:809
Halo(const Pattern &pattern, const int width, const ArrayTypes &... arrays)
Constructor.
Definition Cabana_Grid_Halo.hpp:216
void gather(const ExecutionSpace &exec_space, const ArrayTypes &... arrays) const
Gather data into our ghosts from their owners.
Definition Cabana_Grid_Halo.hpp:288
static KOKKOS_INLINE_FUNCTION void unpackOp(ScatterReduce::Replace, const T &buffer_val, T &array_val)
Reduce an element into the buffer. Replace reduction.
Definition Cabana_Grid_Halo.hpp:800
static KOKKOS_INLINE_FUNCTION std::enable_if_t< 4==ArrayView::rank, void > packElement(const Kokkos::View< char *, memory_space > &buffer, const Kokkos::View< int **, memory_space > &steering, const int element_idx, const ArrayView &array_view)
Definition Cabana_Grid_Halo.hpp:680
static KOKKOS_INLINE_FUNCTION void unpackArray(const ReduceOp &reduce_op, const Kokkos::View< char *, memory_space > &buffer, const Kokkos::View< int **, memory_space > &steering, const int element_idx, const std::integral_constant< std::size_t, 0 >, const Cabana::ParameterPack< ArrayViews... > &array_views)
Unpack an array from a buffer.
Definition Cabana_Grid_Halo.hpp:853
static KOKKOS_INLINE_FUNCTION void unpackOp(ScatterReduce::Min, const T &buffer_val, T &array_val)
Reduce an element into the buffer. Min reduction.
Definition Cabana_Grid_Halo.hpp:782
MPI_Comm getComm(const Array_t &array, const ArrayTypes &... arrays) const
Get the communicator and check to make sure all are the same.
Definition Cabana_Grid_Halo.hpp:464
static KOKKOS_INLINE_FUNCTION void unpackOp(ScatterReduce::Max, const T &buffer_val, T &array_val)
Reduce an element into the buffer. Max reduction.
Definition Cabana_Grid_Halo.hpp:791
static KOKKOS_INLINE_FUNCTION std::enable_if_t< 3==ArrayView::rank, void > packElement(const Kokkos::View< char *, memory_space > &buffer, const Kokkos::View< int **, memory_space > &steering, const int element_idx, const ArrayView &array_view)
Definition Cabana_Grid_Halo.hpp:698
void buildCommData(DecompositionTag decomposition_tag, const int width, const std::array< int, NumSpaceDim > &nid, std::vector< Kokkos::View< char *, memory_space > > &buffers, std::vector< Kokkos::View< int **, memory_space > > &steering, const ArrayTypes &... arrays)
Build communication data.
Definition Cabana_Grid_Halo.hpp:510
static KOKKOS_INLINE_FUNCTION void packArray(const Kokkos::View< char *, memory_space > &buffer, const Kokkos::View< int **, memory_space > &steering, const int element_idx, const std::integral_constant< std::size_t, N >, const Cabana::ParameterPack< ArrayViews... > &array_views)
Pack an array into a buffer.
Definition Cabana_Grid_Halo.hpp:730
MPI_Comm getComm(const Array_t &array) const
Get the communicator.
Definition Cabana_Grid_Halo.hpp:457
void unpackBuffer(const ReduceOp &reduce_op, const ExecutionSpace &exec_space, const Kokkos::View< char *, memory_space > &buffer, const Kokkos::View< int **, memory_space > &steering, ArrayViews... array_views) const
Unpack arrays from a buffer.
Definition Cabana_Grid_Halo.hpp:892
void packBuffer(const ExecutionSpace &exec_space, const Kokkos::View< char *, memory_space > &buffer, const Kokkos::View< int **, memory_space > &steering, ArrayViews... array_views) const
Pack arrays into a buffer.
Definition Cabana_Grid_Halo.hpp:751
auto getLocalGrid(const Array_t &array, const ArrayTypes &... arrays)
Definition Cabana_Grid_Halo.hpp:490
static KOKKOS_INLINE_FUNCTION std::enable_if_t< 3==ArrayView::rank, void > unpackElement(const ReduceOp &reduce_op, const Kokkos::View< char *, memory_space > &buffer, const Kokkos::View< int **, memory_space > &steering, const int element_idx, const ArrayView &array_view)
Definition Cabana_Grid_Halo.hpp:832
void buildSteeringVector(const std::array< IndexSpace< 4 >, NumArray > &spaces, const std::array< std::size_t, NumArray > &value_byte_sizes, const int buffer_bytes, const int buffer_num_element, std::vector< Kokkos::View< int **, memory_space > > &steering)
Build 3d steering vector.
Definition Cabana_Grid_Halo.hpp:555
auto getLocalGrid(const Array_t &array)
Definition Cabana_Grid_Halo.hpp:482
static KOKKOS_INLINE_FUNCTION void unpackArray(const ReduceOp reduce_op, const Kokkos::View< char *, memory_space > &buffer, const Kokkos::View< int **, memory_space > &steering, const int element_idx, const std::integral_constant< std::size_t, N >, const Cabana::ParameterPack< ArrayViews... > &array_views)
Unpack an array from a buffer.
Definition Cabana_Grid_Halo.hpp:869
void buildSteeringVector(const std::array< IndexSpace< 3 >, NumArray > &spaces, const std::array< std::size_t, NumArray > &value_byte_sizes, const int buffer_bytes, const int buffer_num_element, std::vector< Kokkos::View< int **, memory_space > > &steering)
Build 2d steering vector.
Definition Cabana_Grid_Halo.hpp:619
Structured index space.
Definition Cabana_Grid_IndexSpace.hpp:37
Definition Cabana_Grid_Halo.hpp:76
Core: particle data structures and algorithms.
Definition Cabana_AoSoA.hpp:36
auto size(SliceType slice, typename std::enable_if< is_slice< SliceType >::value, int >::type *=0)
Check slice size (differs from Kokkos View).
Definition Cabana_Slice.hpp:1012
ParameterPack< Types... > makeParameterPack(const Types &... ts)
Create a parameter pack.
Definition Cabana_ParameterPack.hpp:189
KOKKOS_FORCEINLINE_FUNCTION std::enable_if< is_parameter_pack< ParameterPack_t >::value, typenameParameterPack_t::templatevalue_type< N > & >::type get(ParameterPack_t &pp)
Get an element from a parameter pack.
Definition Cabana_ParameterPack.hpp:129
Infer array memory space.
Definition Cabana_Grid_Halo.hpp:941
typename ArrayT::memory_space type
Memory space.
Definition Cabana_Grid_Halo.hpp:943
Ghosted decomposition tag.
Definition Cabana_Grid_Types.hpp:197
Owned decomposition tag.
Definition Cabana_Grid_Types.hpp:190
Definition Cabana_Grid_Halo.hpp:173
Definition Cabana_Grid_Halo.hpp:167
Definition Cabana_Grid_Halo.hpp:181
Sum values from neighboring ranks into this rank's data.
Definition Cabana_Grid_Halo.hpp:161
Definition Cabana_ParameterPack.hpp:86